TABLE OF CONTENTS


ABINIT/m_nctk [ Modules ]

[ Top ] [ Modules ]

NAME

 m_nctk

FUNCTION

  Tools and wrappers for NETCDF-IO.

COPYRIGHT

  Copyright (C) 2008-2024 ABINIT group (MG)
  This file is distributed under the terms of the
  GNU General Public License, see ~abinit/COPYING
  or http://www.gnu.org/copyleft/gpl.txt .

TODO

   Remove create_nc_file, write_var_netcdf, the output of OUT.nc is dangereous
     because we can create too many dimensions and get
    nf90_def_dim - NetCDF library returned:   NetCDF: NC_MAX_DIMS exceeded
   Moreover the multiple calls to redef render the IO very inefficient
   That part should be rationalized!

SOURCE

23 #if defined HAVE_CONFIG_H
24 #include "config.h"
25 #endif
26 
27 #include "abi_common.h"
28 
29 MODULE m_nctk
30 
31  use defs_basis
32  use m_abicore
33  use m_build_info
34  use m_errors
35  use, intrinsic :: iso_c_binding
36  use m_xmpi
37 #ifdef HAVE_NETCDF
38  use netcdf
39 #endif
40 
41  use m_fstrings,  only : itoa, sjoin, lstrip, char_count, strcat, endswith, startswith, ltoa
42  use m_io_tools,  only : pick_aname, delete_file, file_exists
43  use m_yaml,      only : DTSET_IDX
44 
45  implicit none
46 
47  private

ABINIT/write_eig [ Functions ]

[ Top ] [ Functions ]

NAME

 write_eig

FUNCTION

 Write the eigenvalues band by band and k point by k point
 in a NetCDF file format

INPUTS

 filname = Filename of the file where the history will be stored

OUTPUT

  (only writing)

SOURCE

2893 subroutine write_eig(eigen,fermie,filename,kptns,mband,nband,nkpt,nsppol,&
2894 & shiftfactor_extfpmd) ! Optional arguments
2895 
2896 !Arguments ------------------------------------
2897 !scalars
2898  character(len=fnlen),intent(in) :: filename
2899  integer,intent(in) :: nkpt,nsppol,mband
2900  real(dp),intent(in) :: fermie
2901  real(dp),optional,intent(in) :: shiftfactor_extfpmd
2902 !arrays
2903  integer,intent(in) :: nband(nkpt*nsppol)
2904  real(dp),intent(in) :: eigen(mband*nkpt*nsppol)
2905  real(dp),intent(in) :: kptns(3,nkpt)
2906 
2907 !Local variables-------------------------------
2908 !scalars
2909  integer :: ncerr,ncid,ii, cmode
2910  integer :: xyz_id,nkpt_id,mband_id,nsppol_id
2911  integer :: eig_id,fermie_id,kpt_id,nbk_id,nbk
2912  integer :: shiftfactor_extfpmd_id
2913  integer :: ikpt,isppol,nband_k,band_index
2914  real(dp):: convrt
2915 !arrays
2916  integer :: dimEIG(3),dimKPT(2),dimNBK(2)
2917  integer :: count2(2),start2(2)
2918  integer :: count3(3),start3(3)
2919  integer :: dim0(0)
2920  real(dp):: band(mband)
2921 
2922 ! *********************************************************************
2923 
2924 #if defined HAVE_NETCDF
2925 
2926  convrt=1.0_dp
2927 
2928 !1. Create netCDF file
2929  !ncerr = nf90_create(path=trim(filename),cmode=NF90_CLOBBER, ncid=ncid)
2930  cmode = def_cmode_for_seq_create
2931  ncerr = nf90_create(path=trim(filename), cmode=cmode, ncid=ncid)
2932  NCF_CHECK_MSG(ncerr," create netcdf EIG file")
2933 
2934 !2. Define dimensions
2935  ncerr = nf90_def_dim(ncid,"xyz",3,xyz_id)
2936  NCF_CHECK_MSG(ncerr," define dimension xyz")
2937 
2938  ncerr = nf90_def_dim(ncid,"mband",mband,mband_id)
2939  NCF_CHECK_MSG(ncerr," define dimension mband")
2940 
2941  ncerr = nf90_def_dim(ncid,"nkpt",nkpt,nkpt_id)
2942  NCF_CHECK_MSG(ncerr," define dimension nkpt")
2943 
2944  ncerr = nf90_def_dim(ncid,"nsppol",nsppol,nsppol_id)
2945  NCF_CHECK_MSG(ncerr," define dimension nsppol")
2946 
2947 !Dimensions for EIGENVALUES
2948  dimEIG = (/ mband_id, nkpt_id, nsppol_id /)
2949 !Dimensions for kpoint positions
2950  dimKPT = (/ xyz_id, nkpt_id /)
2951 !Dimensions for number kpoints per band and spin
2952 !dimNBK = (/ nkpt_id, nsppol_id /)
2953  dimNBK = (/ nkpt_id, nsppol_id /)
2954 
2955 !3. Define variables
2956  call ab_define_var(ncid,dim0,fermie_id,NF90_DOUBLE,&
2957 & "fermie","Chemical potential","Hartree")
2958  call ab_define_var(ncid, dimEIG, eig_id, NF90_DOUBLE,&
2959 & "Eigenvalues",&
2960 & "Values of eigenvalues",&
2961 & "Hartree")
2962  call ab_define_var(ncid, dimKPT, kpt_id, NF90_DOUBLE,"Kptns",&
2963 & "Positions of K-points in reciprocal space",&
2964 & "Dimensionless")
2965  call ab_define_var(ncid, dimNBK, nbk_id, NF90_INT,"NBandK",&
2966 & "Number of bands per kpoint and Spin",&
2967 & "Dimensionless")
2968  if(present(shiftfactor_extfpmd)) then
2969     call ab_define_var(ncid,dim0,shiftfactor_extfpmd_id,NF90_DOUBLE,&
2970 &    "shiftfactor_extfpmd","Extended FPMD shiftfactor","Hartree")
2971  end if
2972 
2973 !4. End define mode
2974  ncerr = nf90_enddef(ncid)
2975  NCF_CHECK_MSG(ncerr," end define mode")
2976 
2977 !5 Write kpoint positions
2978  do ikpt=1,nkpt
2979    start2 = (/ 1, ikpt /)
2980    count2 = (/ 3, 1 /)
2981    ncerr = nf90_put_var(ncid, kpt_id,&
2982 &   kptns(1:3,ikpt),&
2983 &   start = start2,&
2984 &   count = count2)
2985    NCF_CHECK_MSG(ncerr," write variable kptns")
2986  end do
2987 
2988 !6.1 Write chemical potential
2989  ncerr = nf90_put_var(ncid, fermie_id, fermie)
2990  NCF_CHECK_MSG(ncerr," write variable fermie")
2991 
2992 !6.2 Write extfpmd shiftfactor
2993  if(present(shiftfactor_extfpmd)) then
2994    ncerr = nf90_put_var(ncid, shiftfactor_extfpmd_id, shiftfactor_extfpmd)
2995    NCF_CHECK_MSG(ncerr," write variable shiftfactor_extfpmd")
2996  end if
2997 
2998 !6.3 Write eigenvalues
2999  band_index=0
3000  do isppol=1,nsppol
3001    do ikpt=1,nkpt
3002      nband_k=nband(ikpt+(isppol-1)*nkpt)
3003      start3 = (/ 1, ikpt, isppol /)
3004      count3 = (/ mband, 1, 1 /)
3005      band(:)=zero
3006      do ii=1,nband_k
3007        band(ii)=eigen(band_index+ii)
3008      end do
3009      ncerr = nf90_put_var(ncid, eig_id,&
3010 &     band,&
3011 &     start = start3,&
3012 &     count = count3)
3013      NCF_CHECK_MSG(ncerr," write variable band")
3014 
3015      band_index=band_index+nband_k
3016    end do
3017  end do
3018 
3019 !6.4 Write Number of bands per kpoint and Spin
3020 
3021  do isppol=1,nsppol
3022    do ikpt=1,nkpt
3023      start2 = (/ ikpt, 1 /)
3024      count2 = (/ 1, 1 /)
3025      nbk=nband(ikpt+(isppol-1)*nkpt)
3026      ncerr = nf90_put_var(ncid, nbk_id,&
3027 &     nbk,&
3028 &     start = start2)
3029      NCF_CHECK_MSG(ncerr," write variable nband")
3030    end do
3031  end do
3032 
3033 !7 Close file
3034 
3035  ncerr = nf90_close(ncid)
3036  NCF_CHECK_MSG(ncerr," close netcdf EIG file")
3037 #endif
3038 
3039 end subroutine write_eig

m_nctk/ab_define_var [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

 ab_define_var

FUNCTION

 Write the definition of a variable, including units and mnemonics

INPUTS

 ncid = Identifier of the netcdf dataset
 var_dim_id = Identifier of the Dimensions
 var_id     = Identifier of the variable
 var_mnemo  = String of mnemonics
 var_name   = String with the name of the variable
 var_type   = NetCDF type of variable (NF90_DOUBLE, etc)
 var_units  = String of units

OUTPUT

  (only writing)

SOURCE

1350 subroutine ab_define_var(ncid, var_dim_id, var_id, var_type, var_name, var_mnemo, var_units)
1351 
1352 !Arguments ------------------------------------
1353 !scalars
1354  integer, intent(in) :: ncid
1355  integer, intent(out) :: var_id
1356  character(len=*), intent(in) :: var_mnemo,var_units,var_name
1357  integer,intent(in) :: var_type
1358 !arrays
1359  integer,intent(in) :: var_dim_id(:)
1360 
1361 !Local variables-------------------------------
1362 !scalars
1363  integer :: ncerr
1364 
1365 ! *************************************************************************
1366 
1367 #ifdef HAVE_NETCDF
1368  ncerr = nf90_def_var(ncid, trim(var_name), var_type, var_dim_id, var_id)
1369  NCF_CHECK_MSG(ncerr," define variable "//trim(var_name))
1370 
1371  ncerr = nf90_put_att(ncid, var_id,  "units",trim(var_units))
1372  NCF_CHECK_MSG(ncerr," define attribute for "//trim(var_name))
1373 
1374  ncerr = nf90_put_att(ncid, var_id,  "mnemonics", trim(var_mnemo))
1375  NCF_CHECK_MSG(ncerr," define attribute for "//trim(var_name))
1376 #endif
1377 
1378 end subroutine ab_define_var

m_nctk/bail_if_ncerr [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  bail_if_ncerr

FUNCTION

INPUTS

  ncerr=Netcdf error
  line=line number of the file where problem occurred
  file=name of the f90 file containing the caller

SOURCE

692 logical function bail_if_ncerr(ncerr, file, line) result(bail)
693 
694 !Arguments ------------------------------------
695  integer,intent(in) :: ncerr
696  character(len=*),optional,intent(in) :: file
697  integer,optional,intent(in) :: line
698 
699 ! *********************************************************************
700 
701  bail = (ncerr /= nf90_noerr)
702 
703  if (bail) then
704    einfo%ncerr = ncerr
705    einfo%file = "Subroutine Unknown"; if (present(file)) einfo%file = file
706    einfo%line = 0; if (present(line)) einfo%line = line
707    ! Append Netcdf string to user-defined message.
708    write(einfo%msg,'(2a)')'NetCDF library raised: ',trim(nf90_strerror(ncerr))
709  end if
710 
711 end function bail_if_ncerr

m_nctk/collect_datar [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  collect_datar

FUNCTION

 Collect a real-space MPI-FFT distributed array on each proc.

INPUTS

  ngfft(18)=contain all needed information about 3D FFT (see NOTES at beginning of scfcv)
  cplex=1 if real array, 2 for complex
  nfft=Number of FFT points treated by this MPI proc
  nspden=Second dimension of rhor
  rhor(cplex*nfft,nspden)=Array in real space (MPI-FFT distributed)
  fftn3_distrib(n3)=rank of the processors which own fft planes in 3rd dimension.
  fftn3_local(n3)=local i3 indices
  comm_fft=MPI-FFT communicator
  [master]=MPI rank, Optional. If present, the global array is available only on master node.

OUTPUT

   rhor_glob(cplex*nfft_tot,nspden)=Global array

SOURCE

2349 subroutine collect_datar(ngfft,cplex,nfft,nspden,rhor,comm_fft,fftn3_distrib,ffti3_local,rhor_glob,master)
2350 
2351 !Arguments ------------------------------------
2352 !scalars
2353  integer,intent(in) :: cplex,nfft,nspden,comm_fft
2354  integer,optional,intent(in) :: master
2355 !arrays
2356  integer,intent(in) :: ngfft(18)
2357  integer,intent(in) :: fftn3_distrib(ngfft(3)),ffti3_local(ngfft(3))
2358  real(dp),intent(in) :: rhor(cplex*nfft,nspden)
2359  real(dp),intent(out) :: rhor_glob(cplex*product(ngfft(1:3)),nspden)
2360 
2361 !Local variables-------------------------------
2362  integer :: ispden,i1,i2,i3,me_fft,i3_local,my_fftbase,glob_fftbase
2363  integer :: n1,n2,n3,ierr,nfft_tot
2364 
2365 ! *************************************************************************
2366 
2367  nfft_tot = product(ngfft(1:3)); me_fft = xmpi_comm_rank(comm_fft)
2368  n1 = ngfft(1); n2 = ngfft(2); n3 = ngfft(3)
2369 
2370  if (nfft_tot == nfft) then
2371    ! full rhor on each node, just do a copy
2372    rhor_glob = rhor
2373  else
2374    ! if MPI-FFT we have to gather the global array on each node.
2375    rhor_glob = zero
2376    do ispden=1,nspden
2377      do i3=1,n3
2378        if (me_fft == fftn3_distrib(i3)) then
2379          i3_local = ffti3_local(i3)
2380          do i2=1,n2
2381            my_fftbase =   cplex * ( (i2-1)*n1 + (i3_local-1)*n1*n2 )
2382            glob_fftbase = cplex * ( (i2-1)*n1 + (i3-1)*n1*n2 )
2383            do i1=1,cplex * n1
2384              rhor_glob(i1+glob_fftbase,ispden) = rhor(i1+my_fftbase,ispden)
2385            end do
2386          end do
2387        end if
2388      end do
2389    end do
2390    if (present(master)) then
2391      call xmpi_sum_master(rhor_glob,master,comm_fft,ierr)
2392    else
2393      call xmpi_sum(rhor_glob,comm_fft,ierr)
2394    end if
2395  end if
2396 
2397 end subroutine collect_datar

m_nctk/create_nc_file [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

 create_nc_file

FUNCTION

 Create an NetCDF file including a dimension one definition

INPUTS

OUTPUT

 TODO:
  Remove

SOURCE

2761 subroutine create_nc_file (filename,ncid)
2762 
2763 !Arguments ------------------------------------
2764 !scalars
2765  integer,intent(out) :: ncid
2766 !arrays
2767 character(len=*),intent(in) :: filename
2768 
2769 !Local variables-------------------------------
2770 #if defined HAVE_NETCDF
2771 integer :: one_id
2772 integer :: ncerr, cmode
2773 #endif
2774 
2775 ! *************************************************************************
2776 
2777  ncid = 0
2778 #if defined HAVE_NETCDF
2779  ! Create the NetCDF file
2780  !ncerr = nf90_create(path=filename,cmode=NF90_CLOBBER,ncid=ncid)
2781  cmode = def_cmode_for_seq_create
2782  ncerr = nf90_create(path=filename, cmode=cmode, ncid=ncid)
2783  NCF_CHECK_MSG(ncerr, sjoin('Error while creating:', filename))
2784  ncerr=nf90_def_dim(ncid,'one',1,one_id)
2785  NCF_CHECK_MSG(ncerr,'nf90_def_dim')
2786 #endif
2787 
2788  end subroutine create_nc_file

m_nctk/distrib_datar [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  distrib_datar

FUNCTION

 distribute a real-space MPI-FFT

INPUTS

  ngfft(18)=contain all needed information about 3D FFT (see NOTES at beginning of scfcv)
  cplex=1 if real array, 2 for complex
  nfft=Number of FFT points treated by this MPI proc
  nspden=Second dimension of rhor
  rhor_glob(cplex*nfft_tot,nspden)=Global array
  master=The rank of the node that owns the global array.
  comm_fft=MPI-FFT communicator
  fftn3_distrib(n3)=rank of the processors which own fft planes in 3rd dimension.
  fftn3_local(n3)=local i3 indices

OUTPUT

  rhor(cplex*nfft,nspden)=Array in real space (MPI-FFT distributed)

SOURCE

2425 subroutine distrib_datar(ngfft,cplex,nfft,nspden,rhor_glob,master,comm_fft,fftn3_distrib,ffti3_local,rhor)
2426 
2427 !Arguments ------------------------------------
2428 !scalars
2429  integer,intent(in) :: cplex,nfft,nspden,comm_fft,master
2430 !arrays
2431  integer,intent(in) :: ngfft(18)
2432  integer,intent(in) :: fftn3_distrib(ngfft(3)),ffti3_local(ngfft(3))
2433  real(dp),intent(out) :: rhor(cplex*nfft,nspden)
2434  real(dp),intent(inout) :: rhor_glob(cplex*product(ngfft(1:3)),nspden)
2435 
2436 !Local variables-------------------------------
2437  integer :: ispden,i1,i2,i3,me_fft,i3_local,my_fftbase,glob_fftbase
2438  integer :: n1,n2,n3,ierr,nfft_tot
2439 
2440 ! *************************************************************************
2441 
2442  nfft_tot = product(ngfft(1:3)); me_fft = xmpi_comm_rank(comm_fft)
2443  n1 = ngfft(1); n2 = ngfft(2); n3 = ngfft(3)
2444 
2445  if (nfft_tot == nfft) then
2446    ! full rhor on each node, just do a copy
2447    rhor = rhor_glob
2448  else
2449    ! if MPI-FFT we have to gather the global array on each node.
2450    call xmpi_bcast(rhor_glob,master,comm_fft,ierr)
2451    do ispden=1,nspden
2452      do i3=1,n3
2453        if (me_fft == fftn3_distrib(i3)) then
2454          i3_local = ffti3_local(i3)
2455          do i2=1,n2
2456            my_fftbase =   cplex * ( (i2-1)*n1 + (i3_local-1)*n1*n2 )
2457            glob_fftbase = cplex * ( (i2-1)*n1 + (i3-1)*n1*n2 )
2458            do i1=1,cplex * n1
2459              rhor(i1+my_fftbase,ispden) = rhor_glob(i1+glob_fftbase,ispden)
2460            end do
2461          end do
2462        end if
2463      end do
2464    end do
2465  end if
2466 
2467 end subroutine distrib_datar

m_nctk/ncfdim_t [ Types ]

[ Top ] [ m_nctk ] [ Types ]

NAME

 nctkdim_t

FUNCTION

  Stores the name and the value of a netcdf dimension

SOURCE

119  type,public :: nctkdim_t
120    character(len=nctk_slen) :: name   ! name of the dimension.
121    integer :: value                   ! value of the dimension.
122    !integer :: id
123  end type nctkdim_t

m_nctk/nctk_add_etsf_header [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_add_etsf_header

FUNCTION

   Add the etsf-io header to a file associated to a netcdf file handler.

INPUTS

  ncid=Netcdf file identifier.
  * version = the number of version to be created.
  * title = (optional) a title for the file (80 characters max).
  * history = (optional) the first line of history (1024 characters max).
  * with_etsf_header = (optional) if true, will create a header
                       as defined in the ETSF specifications (default is .true.).
                       When value is .false., arguments title, history and version
                       are ignored.

SOURCE

945 integer function nctk_add_etsf_header(ncid, title, history) result(ncerr)
946 
947 !Arguments ------------------------------------
948  integer,intent(in) :: ncid
949  character(len=*),optional,intent(in) :: title,history
950 
951 !Local variables-------------------------------
952  !integer :: ncerr
953  character(len=*), parameter :: file_format = "ETSF Nanoquanta"
954  character(len=*), parameter :: conventions = "http://www.etsf.eu/fileformats/"
955  real :: format_version = 3.3 ! Real is not a good choice for a version!
956 
957 ! *********************************************************************
958 
959  ncerr = nctk_set_defmode(ncid)
960  if (ncerr /= nf90_noerr) return
961 
962  ! The file format
963  ncerr = nf90_put_att(ncid, NF90_GLOBAL, "file_format", file_format)
964  if (ncerr /= nf90_noerr) return
965 
966  ! The version
967  ncerr = nf90_put_att(ncid, NF90_GLOBAL, "file_format_version", format_version)
968  if (ncerr /= nf90_noerr) return
969 
970  ! The conventions
971  ncerr = nf90_put_att(ncid, NF90_GLOBAL, "Conventions", conventions)
972  if (ncerr /= nf90_noerr) return
973 
974  ! The history if present
975  if (present(history)) then
976    ncerr = nf90_put_att(ncid, NF90_GLOBAL, "history", history(1:min(1024, len(history))))
977    if (ncerr /= nf90_noerr) return
978  end if
979 
980  ! The title if present
981  if (present(title)) then
982    ncerr = nf90_put_att(ncid, NF90_GLOBAL, "title", title(1:min(80, len(title))))
983    if (ncerr /= nf90_noerr) return
984  end if
985 
986  ! These are extensions not in the standard.
987  ! Add info on the code that produced this file
988  ncerr = nf90_put_att(ncid, NF90_GLOBAL, "code", "Abinit")
989  if (ncerr /= nf90_noerr) return
990 
991  ncerr = nf90_put_att(ncid, NF90_GLOBAL, "code_version", ABINIT_VERSION)
992  if (ncerr /= nf90_noerr) return
993 
994 end function nctk_add_etsf_header

m_nctk/nctk_def_array_list [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_def_array_list

FUNCTION

  Define list of arrays with a given type, Return immediately if error

INPUTS

  ncid=Netcdf identifier.
  nctk_arrays(:)=List of array descriptors.
  [defmode]=If True, the nc file is set in define mode (default=False)
  [prefix]=Prefix added to varnames and dimensions. Empty string if not specified.

SOURCE

1666 integer function nctk_def_array_list(ncid, nctk_arrays, defmode, prefix) result(ncerr)
1667 
1668 !Arguments ------------------------------------
1669 !scalars
1670  integer,intent(in) :: ncid
1671  logical,optional,intent(in) :: defmode
1672  character(len=*),optional,intent(in) :: prefix
1673 !arrays
1674  type(nctkarr_t),intent(in) :: nctk_arrays(:)
1675 
1676 !Local variables-------------------------------
1677 !scalars
1678  integer :: ia
1679 
1680 ! *********************************************************************
1681 
1682  ncerr = nf90_noerr
1683  if (present(defmode)) then
1684    if (defmode) then
1685      NCF_CHECK(nctk_set_defmode(ncid))
1686    end if
1687  end if
1688 
1689  do ia=1,size(nctk_arrays)
1690    if (present(prefix)) then
1691      NCF_CHECK(nctk_def_one_array(ncid, nctk_arrays(ia), prefix=prefix))
1692    else
1693      NCF_CHECK(nctk_def_one_array(ncid, nctk_arrays(ia)))
1694    end if
1695  end do
1696 
1697 end function nctk_def_array_list

m_nctk/nctk_def_basedims [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_def_basedims

FUNCTION

  Define the basic dimensions used in ETSF-IO files.

INPUTS

  ncid=Netcdf identifier.
  [defmode]=If True, the nc file is set in define mode (default=False)

SOURCE

1292 integer function nctk_def_basedims(ncid, defmode) result(ncerr)
1293 
1294 !Arguments ------------------------------------
1295  integer,intent(in) :: ncid
1296  logical,optional,intent(in) :: defmode
1297 
1298 ! *********************************************************************
1299 
1300  ncerr = nf90_noerr
1301 
1302  if (present(defmode)) then
1303    if (defmode) then
1304      NCF_CHECK(nctk_set_defmode(ncid))
1305    end if
1306  end if
1307 
1308  ! Basic ETSF-IO dimensions that should be always present in the file.
1309  ncerr = nctk_def_dims(ncid, [&
1310    nctkdim_t("complex", 2), nctkdim_t("symbol_length", 2), nctkdim_t("character_string_length", etsfio_charlen),&
1311    nctkdim_t("number_of_cartesian_directions", 3), nctkdim_t("number_of_reduced_dimensions", 3),&
1312    nctkdim_t("number_of_vectors", 3) &
1313  ])
1314  NCF_CHECK(ncerr)
1315 
1316  ! Useful integers.
1317  ncerr = nctk_def_dims(ncid, [ &
1318    nctkdim_t("one", 1), nctkdim_t("two", 2), nctkdim_t("three", 3), &
1319    nctkdim_t("four", 4), nctkdim_t("five", 5), nctkdim_t("six", 6), &
1320    nctkdim_t("seven", 7), nctkdim_t("eight", 8), nctkdim_t("nine", 9), nctkdim_t("ten", 10), &
1321    nctkdim_t("fnlen", fnlen + 1) &
1322  ])
1323  NCF_CHECK(ncerr)
1324 
1325 end function nctk_def_basedims

m_nctk/nctk_def_dim_list [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_def_dim_list

FUNCTION

  Define list of dimensions variables and write their values.
  Return immediately if error

INPUTS

  ncid=Netcdf identifier.
  dimnames(:)=List of strings with the name of the dimensions
  values(:)=List of integer scalars
  [defmode]=If True, the nc file is set in define mode (default=False)
  [prefix]=Prefix added to varnames and dimensions. Empty string if not specified.

SOURCE

1206 integer function nctk_def_dim_list(ncid, nctkdims, defmode, prefix) result(ncerr)
1207 
1208 !Arguments ------------------------------------
1209 !scalars
1210  integer,intent(in) :: ncid
1211  logical,optional,intent(in) :: defmode
1212  character(len=*),optional,intent(in) :: prefix
1213 !arrays
1214  type(nctkdim_t),intent(in) :: nctkdims(:)
1215 
1216 !Local variables-------------------------------
1217 !scalars
1218  integer :: ii
1219 
1220 ! *********************************************************************
1221 
1222  ncerr = nf90_noerr
1223  if (present(defmode)) then
1224    if (defmode) then
1225      NCF_CHECK(nctk_set_defmode(ncid))
1226    end if
1227  end if
1228 
1229  do ii=1,size(nctkdims)
1230    if (present(prefix)) then
1231      ncerr = nctk_def_one_dim(ncid, nctkdims(ii), prefix=prefix)
1232    else
1233      ncerr = nctk_def_one_dim(ncid, nctkdims(ii))
1234    end if
1235    if (ncerr /= nf90_noerr) return
1236  end do
1237 
1238 end function nctk_def_dim_list

m_nctk/nctk_def_dpscalars [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_def_dpscalars

FUNCTION

  Define list of double precision **scalar** variables. Return immediately if error

INPUTS

  ncid=Netcdf identifier.
  varnames(:)=List of strings with the name of the variables
  [defmode]=If True, the nc file is set in define mode (default=False)
  [prefix]=Prefix added to varnames and dimensions. Empty string if not specified.

SOURCE

1502 integer function nctk_def_dpscalars(ncid, varnames, defmode, prefix) result(ncerr)
1503 
1504 !Arguments ------------------------------------
1505 !scalars
1506  integer,intent(in) :: ncid
1507  logical,optional,intent(in) :: defmode
1508  character(len=*),optional,intent(in) :: prefix
1509 !arrays
1510  character(len=*),intent(in) :: varnames(:)
1511 
1512 !Local variables-------------------------------
1513  character(len=nctk_slen) :: prefix_
1514 
1515 ! *********************************************************************
1516  prefix_ = ""; if (present(prefix)) prefix_ = prefix
1517 
1518  if (present(defmode)) then
1519    ncerr = nctk_def_scalars_type(ncid, varnames, nf90_double, defmode=defmode, prefix=prefix_)
1520  else
1521    ncerr = nctk_def_scalars_type(ncid, varnames, nf90_double, prefix=prefix_)
1522  end if
1523 
1524 end function nctk_def_dpscalars

m_nctk/nctk_def_iscalars [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_def_iscalars

FUNCTION

  Define list of integer **scalar** variables. Return immediately if error

INPUTS

  ncid=Netcdf identifier.
  varnames(:)=List of strings with the name of the variables
  [defmode]=If True, the nc file is set in define mode (default=False)
  [prefix]=Prefix added to varnames and dimensions. Empty string if not specified.

SOURCE

1462 integer function nctk_def_iscalars(ncid, varnames, defmode, prefix) result(ncerr)
1463 
1464 !Arguments ------------------------------------
1465 !scalars
1466  integer,intent(in) :: ncid
1467  logical,optional,intent(in) :: defmode
1468  character(len=*),optional,intent(in) :: prefix
1469 !arrays
1470  character(len=*),intent(in) :: varnames(:)
1471 
1472 ! *********************************************************************
1473 
1474  if (present(defmode)) then
1475    ncerr = nctk_def_scalars_type(ncid, varnames, nf90_int, defmode=defmode)
1476  else
1477    if (present(prefix)) then
1478      ncerr = nctk_def_scalars_type(ncid, varnames, nf90_int, prefix=prefix)
1479    else
1480      ncerr = nctk_def_scalars_type(ncid, varnames, nf90_int)
1481    end if
1482  end if
1483 
1484 end function nctk_def_iscalars

m_nctk/nctk_def_one_array [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_def_one_array

FUNCTION

  Define list of arrays with a given type, Return immediately if error

INPUTS

  ncid=Netcdf identifier.
  nctk_array=Array descriptor.
  [defmode]=If True, the nc file is set in define mode (default=False)
  [prefix]=Prefix added to varnames and dimensions. Empty string if not specified.

SOURCE

1543 integer function nctk_def_one_array(ncid, nctk_array, defmode, varid, prefix) result(ncerr)
1544 
1545 !Arguments ------------------------------------
1546 !scalars
1547  integer,intent(in) :: ncid
1548  logical,optional,intent(in) :: defmode
1549  integer,optional,intent(out) :: varid
1550  type(nctkarr_t),intent(in) :: nctk_array
1551  character(len=*),optional,intent(in) :: prefix
1552 
1553 !Local variables-------------------------------
1554 !scalars
1555  integer :: ii,xtype,prev,cnt,nn,vid
1556  character(len=500) :: msg
1557 !arrays
1558  integer :: dimids(NF90_MAX_DIMS),dimvals(NF90_MAX_DIMS)
1559  character(len=nctk_slen) :: sarr(NF90_MAX_DIMS), string, pre, vname, dimname
1560  type(nctkvar_t) :: var
1561 
1562 ! *********************************************************************
1563 
1564  pre = ""; if (present(prefix)) pre = prefix
1565 
1566  ncerr = nf90_noerr
1567  if (present(defmode)) then
1568    if (defmode) then
1569      NCF_CHECK(nctk_set_defmode(ncid))
1570    end if
1571  end if
1572 
1573  xtype = str2xtype(nctk_array%dtype)
1574  vname = strcat(pre, nctk_array%name)
1575 
1576  ! Build array of strings with the dimensions.
1577  string = nctk_array%shape_str
1578  nn = char_count(string, ",")
1579  ABI_CHECK(nn <= NF90_MAX_DIMS, "Too many dimensions!")
1580 
1581  ! Parse dimension names and add prefix (if any).
1582  if (nn == 0) then
1583    cnt = 1
1584    dimname = lstrip(string)
1585    if (any(dimname == NCTK_IMPLICIT_DIMS)) then
1586      sarr(1) = dimname
1587    else
1588      sarr(1) = strcat(pre, dimname)
1589    end if
1590  else
1591    prev = 0; cnt = 0
1592    do ii=1,len_trim(string)
1593      if (string(ii:ii) == ",") then
1594        cnt = cnt + 1
1595        dimname = lstrip(string(prev+1:ii-1))
1596        if (any(dimname == NCTK_IMPLICIT_DIMS)) then
1597          sarr(cnt) = dimname
1598        else
1599          sarr(cnt) = strcat(pre, dimname)
1600        end if
1601        prev = ii
1602      end if
1603    end do
1604    cnt = cnt + 1
1605    dimname = lstrip(string(prev+1:ii-1))
1606    if (any(dimname == NCTK_IMPLICIT_DIMS)) then
1607      sarr(cnt) = dimname
1608    else
1609      sarr(cnt) = strcat(pre, dimname)
1610    end if
1611  end if
1612 
1613  ! Get dimids
1614  do ii=1,cnt
1615    NCF_CHECK_MSG(nf90_inq_dimid(ncid, sarr(ii), dimids(ii)), sarr(ii))
1616    NCF_CHECK(nf90_inquire_dimension(ncid, dimids(ii), len=dimvals(ii)))
1617  end do
1618 
1619  ! Check if dimension already exists.
1620  ! Variable already exists. Check if type and dimensions agree
1621  if (nf90_inq_varid(ncid, vname, vid) == nf90_noerr)  then
1622    call var_from_id(ncid, vid, var)
1623    if (.not. (var%xtype == xtype .and. var%ndims == cnt)) then
1624       write(msg,"(4a,2(2(a,i0),a))")&
1625         "variable ",trim(vname)," already exists with a different definition:",ch10,&
1626         "In file:     xtype = ",var%xtype,", ndims = ",var%ndims,ch10,&
1627         "From caller: xtype = ",xtype,", ndims = ",cnt,ch10
1628       ABI_ERROR(msg)
1629    end if
1630    if (any(dimvals(1:cnt) /= var%dimlens(1:var%ndims))) then
1631       write(msg,"(4a,2(3a))")&
1632         "variable ",trim(vname)," already exists but with different shape.",ch10,&
1633         "In file:     dims = ",trim(ltoa(var%dimlens(:var%ndims))),ch10,&
1634         "From caller  dims = ",trim(ltoa(dimvals(:cnt))),ch10
1635       ABI_ERROR(msg)
1636    end if
1637    if (present(varid)) varid = vid
1638    return
1639  end if
1640 
1641  ! Define variable since it doesn't exist.
1642  ncerr = nf90_def_var(ncid, vname, xtype, dimids(1:cnt), vid)
1643  NCF_CHECK(ncerr)
1644 
1645  if (present(varid)) varid = vid
1646 
1647 end function nctk_def_one_array

m_nctk/nctk_def_one_dim [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_def_one_dim

FUNCTION

  Define list of dimensions variables and write their values.
  Return immediately if error

INPUTS

  ncid=Netcdf identifier.
  dimnames(:)=List of strings with the name of the dimensions
  values(:)=List of integer scalars
  [defmode]=If True, the nc file is set in define mode (default=False)
  [prefix]=Prefix added to varnames and dimensions. Empty string if not specified.

SOURCE

1134 integer function nctk_def_one_dim(ncid, nctkdim, defmode, prefix) result(ncerr)
1135 
1136 !Arguments ------------------------------------
1137 !scalars
1138  integer,intent(in) :: ncid
1139  logical,optional,intent(in) :: defmode
1140  character(len=*),optional,intent(in) :: prefix
1141 !arrays
1142  type(nctkdim_t),intent(in) :: nctkdim
1143 
1144 !Local variables-------------------------------
1145  integer :: dimid,dimlen
1146  character(len=nctk_slen) :: dname
1147  character(len=500) :: msg
1148 ! *********************************************************************
1149 
1150  ncerr = nf90_noerr
1151 
1152  if (present(defmode)) then
1153    if (defmode) then
1154      NCF_CHECK(nctk_set_defmode(ncid))
1155    end if
1156  end if
1157 
1158  if (present(prefix)) then
1159    if (any(nctkdim%name == NCTK_IMPLICIT_DIMS)) then
1160      dname = nctkdim%name
1161    else
1162      dname = strcat(prefix, nctkdim%name)
1163    end if
1164  else
1165    dname = nctkdim%name
1166  end if
1167 
1168  ! if dimension already exists, test whether it has the same value else define new dim.
1169  ncerr = nf90_inq_dimid(ncid, dname, dimid)
1170 
1171  if (ncerr == nf90_noerr) then
1172    NCF_CHECK(nf90_inquire_dimension(ncid, dimid, len=dimlen))
1173    if (dimlen /= nctkdim%value) then
1174      write(msg, "(4a,2(a,i0))")&
1175         "dimension ", trim(dname)," already exists but with a different value",ch10,&
1176         "from file: ", dimlen, "; about to write: ", nctkdim%value
1177      ABI_ERROR(msg)
1178    end if
1179  else
1180    ncerr = nf90_def_dim(ncid, dname, nctkdim%value, dimid)
1181    NCF_CHECK(ncerr)
1182  end if
1183 
1184 end function nctk_def_one_dim

m_nctk/nctk_def_scalars_type [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_def_scalars_type

FUNCTION

  Define list of **scalar** variables with a given type, Return immediately if error

INPUTS

  ncid=Netcdf identifier.
  varnames(:)=List of strings with the name of the variables
  xtype=Type of the variables
  [defmode]=If True, the nc file is set in define mode (default=False)
  [prefix]=Prefix added to varnames and dimensions. Empty string if not specified.

SOURCE

1397 integer function nctk_def_scalars_type(ncid, varnames, xtype, defmode, prefix) result(ncerr)
1398 
1399 !Arguments ------------------------------------
1400 !scalars
1401  integer,intent(in) :: ncid,xtype
1402  logical,optional,intent(in) :: defmode
1403  character(len=*),optional,intent(in) :: prefix
1404 !arrays
1405  character(len=*),intent(in) :: varnames(:)
1406 
1407 !Local variables-------------------------------
1408 !scalars
1409  integer :: ii,varid
1410  character(len=nctk_slen) :: vname
1411  type(nctkvar_t) :: var
1412 
1413 ! *********************************************************************
1414 
1415  ncerr = nf90_noerr
1416  if (present(defmode)) then
1417    if (defmode) then
1418      NCF_CHECK(nctk_set_defmode(ncid))
1419    end if
1420  end if
1421 
1422  ! Special case where dimension is null
1423  do ii=1,size(varnames)
1424    vname = varnames(ii)
1425    if (present(prefix)) vname = strcat(prefix, vname)
1426 
1427    if (nf90_inq_varid(ncid, vname, varid) == nf90_noerr)  then
1428        ! Variable already exists. Check if type and dimensions agree
1429        call var_from_id(ncid, varid, var)
1430 
1431        if (.not. (var%xtype == xtype .and. var%ndims == 0)) then
1432          ABI_ERROR("variable already exists with a different definition.")
1433        else
1434          cycle ! Dimension matches, skip definition.
1435        end if
1436 
1437    else
1438      ! Define variable since it doesn't exist.
1439      ncerr = nf90_def_var(ncid, vname, xtype, varid)
1440      NCF_CHECK(ncerr)
1441    end if
1442  end do
1443 
1444 end function nctk_def_scalars_type

m_nctk/nctk_defnwrite_dpvars [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_defnwrite_dpvars

FUNCTION

  Define list of real(dp) **scalar** variables and write their values.
  Return immediately if error

INPUTS

  ncid=Netcdf identifier.
  varnames(:)=List of strings with the name of the variables
  values(:)=List of integer scalars

SOURCE

1859 integer function nctk_defnwrite_dpvars(ncid, varnames, values) result(ncerr)
1860 
1861 !Arguments ------------------------------------
1862 !scalars
1863  integer,intent(in) :: ncid
1864 !arrays
1865  real(dp),intent(in) :: values(:)
1866  character(len=*),intent(in) :: varnames(:)
1867 
1868 !Local variables-------------------------------
1869 !scalars
1870  integer :: ii,varid
1871 !arrays
1872 
1873 ! *********************************************************************
1874  ncerr = nf90_noerr
1875 
1876  ABI_CHECK(size(varnames) == size(values), "Different size in varnames, values")
1877 
1878  ncerr = nctk_def_dpscalars(ncid, varnames, defmode=.True.)
1879  NCF_CHECK(ncerr)
1880 
1881  NCF_CHECK(nctk_set_datamode(ncid))
1882  do ii=1,size(varnames)
1883    !write(std_out,*)varnames(ii)
1884    varid = nctk_idname(ncid, varnames(ii))
1885    NCF_CHECK(nf90_put_var(ncid, varid, values(ii)))
1886  end do
1887 
1888 end function nctk_defnwrite_dpvars

m_nctk/nctk_defnwrite_ivars [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_defnwrite_ivars

FUNCTION

  Define list of integer **scalar** variables and write their values.
  Return immediately if error

INPUTS

  ncid=Netcdf identifier.
  varnames(:)=List of strings with the name of the variables
  values(:)=List of integer scalars

SOURCE

1815 integer function nctk_defnwrite_ivars(ncid, varnames, values) result(ncerr)
1816 
1817 !Arguments ------------------------------------
1818 !scalars
1819  integer,intent(in) :: ncid
1820 !arrays
1821  integer,intent(in) :: values(:)
1822  character(len=*),intent(in) :: varnames(:)
1823 
1824 !Local variables-------------------------------
1825 !scalars
1826  integer :: ii,varid
1827 
1828 ! *********************************************************************
1829 
1830  ABI_CHECK(size(varnames) == size(values), "Different size in varnames, values")
1831 
1832  ncerr = nctk_def_iscalars(ncid, varnames, defmode=.True.)
1833  NCF_CHECK(ncerr)
1834 
1835  NCF_CHECK(nctk_set_datamode(ncid))
1836  do ii=1,size(varnames)
1837    varid = nctk_idname(ncid, varnames(ii))
1838    NCF_CHECK(nf90_put_var(ncid, varid, values(ii)))
1839  end do
1840 
1841 end function nctk_defnwrite_ivars

m_nctk/nctk_defwrite_nonana_raman_terms [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

 nctk_defwrite_nonana_raman_terms

FUNCTION

 Write the Raman susceptiblities for q-->0 along different directions in the netcdf file.

INPUTS

  ncid=netcdf file id.
  iphl2=Index of the q-point to be written to file.
  nph2l=Number of qpoints.
  rsus(3*natom,3,3)=List of Raman susceptibilities along the direction corresponding to iphl2.
  natom=Number of atoms

OUTPUT

  Only writing.

SOURCE

2653 subroutine nctk_defwrite_nonana_raman_terms(ncid, iphl2, nph2l, natom, rsus, mode)
2654 
2655 !Arguments ------------------------------------
2656 !scalars
2657  integer,intent(in) :: ncid,natom,iphl2,nph2l
2658  character(len=*),intent(in) :: mode
2659 !arrays
2660  real(dp),intent(in) :: rsus(3*natom,3,3)
2661 
2662 !Local variables-------------------------------
2663 !scalars
2664  integer :: ncerr, raman_sus_varid
2665 
2666 ! *************************************************************************
2667 
2668 !Fake use of nph2l, to keep it as argument. This should be removed when nph2l will be used.
2669  if(.false.)then
2670   ncerr=nph2l
2671  endif
2672 
2673  select case (mode)
2674  case ("define")
2675    NCF_CHECK(nctk_def_basedims(ncid, defmode=.True.))
2676    ncerr = nctk_def_arrays(ncid, [ nctkarr_t("non_analytical_raman_sus", "dp", &
2677 "number_of_non_analytical_directions,number_of_phonon_modes,number_of_cartesian_directions,number_of_cartesian_directions")])
2678    NCF_CHECK(ncerr)
2679 
2680    NCF_CHECK(nctk_set_datamode(ncid))
2681 
2682  case ("write")
2683    NCF_CHECK(nf90_inq_varid(ncid, "non_analytical_raman_sus", raman_sus_varid))
2684    ncerr = nf90_put_var(ncid,raman_sus_varid,rsus,&
2685      start=[iphl2,1,1,1], count=[1,3*natom,3,3])
2686    NCF_CHECK(ncerr)
2687 
2688  case default
2689    ABI_ERROR(sjoin("Wrong value for mode", mode))
2690  end select
2691 
2692 end subroutine nctk_defwrite_nonana_raman_terms

m_nctk/nctk_defwrite_nonana_terms [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

 nctk_defwrite_nonana_terms

FUNCTION

  Write to ncfile the phonon frequencies and displacements for q --> 0 in the presence of non-analytical behaviour.

INPUTS

  ncid=netcdf file id.
  iphl2=Index of the q-point to be written to file
  nph2l=Number of qpoints.
  qph2l(3,nph2l)=List of phonon wavevector directions along which the non-analytical correction
    to the Gamma-point phonon frequencies will be calculated
    The direction is in CARTESIAN COORDINATES
  natom=Number of atoms
  phfrq(3*natom)=Phonon frequencies in Ha
  cart_displ(2,3*natom,3*natom)=displacements in CARTESIAN coordinates.

OUTPUT

  Only writing.

SOURCE

2585 subroutine nctk_defwrite_nonana_terms(ncid, iphl2, nph2l, qph2l, natom, phfrq, cart_displ, mode)
2586 
2587 !Arguments ------------------------------------
2588 !scalars
2589  integer,intent(in) :: ncid,iphl2,nph2l,natom
2590  character(len=*),intent(in) :: mode
2591 !arrays
2592  real(dp),intent(in) :: qph2l(3, nph2l)
2593  real(dp),intent(in) :: phfrq(3*natom)
2594  real(dp),intent(in) :: cart_displ(2,3*natom,3*natom)
2595 
2596 !Local variables-------------------------------
2597 !scalars
2598  integer :: ncerr, na_phmodes_varid, na_phdispl_varid
2599 
2600 ! *************************************************************************
2601 
2602  select case (mode)
2603  case ("define")
2604    !NCF_CHECK(nctk_def_basedims(ncid, defmode=.True.))
2605    ncerr = nctk_def_dims(ncid, [nctkdim_t("number_of_non_analytical_directions", nph2l)], defmode=.True.)
2606    NCF_CHECK(ncerr)
2607 
2608    ncerr = nctk_def_arrays(ncid, [&
2609      nctkarr_t('non_analytical_directions', "dp", "number_of_cartesian_directions, number_of_non_analytical_directions"),&
2610      nctkarr_t('non_analytical_phonon_modes', "dp", "number_of_phonon_modes, number_of_non_analytical_directions"),&
2611      nctkarr_t('non_analytical_phdispl_cart', "dp", &
2612                "two, number_of_phonon_modes, number_of_phonon_modes, number_of_non_analytical_directions")])
2613    NCF_CHECK(ncerr)
2614 
2615    NCF_CHECK(nctk_set_datamode(ncid))
2616    NCF_CHECK(nf90_put_var(ncid, nctk_idname(ncid, "non_analytical_directions"), qph2l))
2617 
2618  case ("write")
2619 
2620    NCF_CHECK(nf90_inq_varid(ncid, "non_analytical_phonon_modes", na_phmodes_varid))
2621    NCF_CHECK(nf90_put_var(ncid,na_phmodes_varid,phfrq*Ha_eV,start=[1, iphl2], count=[3*natom, 1]))
2622    NCF_CHECK(nf90_inq_varid(ncid, "non_analytical_phdispl_cart", na_phdispl_varid))
2623    ncerr = nf90_put_var(ncid,na_phdispl_varid,cart_displ*Bohr_Ang,&
2624    start=[1,1,1,iphl2], count=[2,3*natom,3*natom, 1])
2625    NCF_CHECK(ncerr)
2626 
2627  case default
2628    ABI_ERROR(sjoin("Wrong value for mode", mode))
2629  end select
2630 
2631 end subroutine nctk_defwrite_nonana_terms

m_nctk/nctk_defwrite_raman_terms [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

 nctk_defwrite_raman_terms

FUNCTION

 Write the Raman susceptiblities for q=0 and also the phonon frequncies at gamma.

INPUTS

  ncid=netcdf file id.
  rsus(3*natom,3,3)=List of Raman susceptibilities.
  natom=Number of atoms

OUTPUT

  Only writing.

SOURCE

2712 subroutine nctk_defwrite_raman_terms(ncid, natom, rsus, phfrq)
2713 
2714 !Arguments ------------------------------------
2715 !scalars
2716  integer,intent(in) :: ncid,natom
2717 !arrays
2718  real(dp),intent(in) :: rsus(3*natom,3,3)
2719  real(dp),intent(in) :: phfrq(3*natom)
2720 
2721 !Local variables-------------------------------
2722 !scalars
2723  integer :: ncerr, raman_sus_varid, phmodes_varid
2724 
2725 ! *************************************************************************
2726 
2727  NCF_CHECK(nctk_def_basedims(ncid, defmode=.True.))
2728  ncerr = nctk_def_arrays(ncid, [ nctkarr_t("raman_sus", "dp", &
2729   "number_of_phonon_modes,number_of_cartesian_directions,number_of_cartesian_directions"), &
2730   nctkarr_t("gamma_phonon_modes", "dp", "number_of_phonon_modes")])
2731  NCF_CHECK(ncerr)
2732 
2733  NCF_CHECK(nctk_set_datamode(ncid))
2734 
2735  NCF_CHECK(nf90_inq_varid(ncid, "raman_sus", raman_sus_varid))
2736  NCF_CHECK(nf90_put_var(ncid,raman_sus_varid,rsus))
2737  NCF_CHECK(nf90_inq_varid(ncid, "gamma_phonon_modes", phmodes_varid))
2738  NCF_CHECK(nf90_put_var(ncid,phmodes_varid,phfrq*Ha_eV))
2739 
2740 end subroutine nctk_defwrite_raman_terms

m_nctk/nctk_fort_or_ncfile [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_fort_or_ncfile

FUNCTION

  Return the iomode used to perform IO operations on filename.
  If filename does not exist, a similar file with extension `.nc` is tried
  and iomode is set to IO_MODE_ETSF if the file exists.
  This trick is used to run the Abinit test suite in netcdf mode without changing the input files.
  The modification (if any) is logged to std_out.

SIDE EFFECTS

  filename=Tentative filename in input. Changed to netcdf file if input filename does not exist
   and a file with netcdf extension is found.

OUTPUT

  iomode=Flag selecting the IO library. Set to IO_MODE_ETSF if netcdf file, else IO_MODE_MPI
    if MPI supports it, finally IO_MODE_FORTRAN
  errmsg=String with error message. Use `if (len_trim(errmsg) /= 0) ABI_ERROR(errmsg)`
    to handle possible errors in the caller.

SOURCE

465 subroutine nctk_fort_or_ncfile(filename, iomode, errmsg)
466 
467  character(len=*),intent(inout) :: filename
468  character(len=*),intent(out) :: errmsg
469  integer,intent(out) :: iomode
470 
471 ! *********************************************************************
472   errmsg = ""
473 
474  ! Default value
475 #ifdef HAVE_MPI_IO
476  iomode = IO_MODE_MPI
477 #else
478  iomode = IO_MODE_FORTRAN
479 #endif
480 
481  !  Checking the existence of data file
482  if (.not.file_exists(filename)) then
483    ! Trick needed to run Abinit test suite in netcdf mode.
484    if (file_exists(nctk_ncify(filename))) then
485      write(std_out,"(3a)")"- File: ",trim(filename)," does not exist but found netcdf file with similar name."
486      filename = nctk_ncify(filename); iomode = IO_MODE_ETSF
487    end if
488    if (.not. file_exists(filename)) then
489      errmsg = 'Missing file: '//trim(filename)
490    end if
491  end if
492 
493 end subroutine nctk_fort_or_ncfile

m_nctk/nctk_get_dim [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_get_dim

FUNCTION

  Get the value of a dimension from its name.

INPUTS

  ncid=Netcdf identifier.
  dimname=Name of the dimension.
  [datamode]=If True, the nc file is set in data mode (default=False)

OUTPUT

  dimlen=Value of the dimension.

SOURCE

1966 integer function nctk_get_dim(ncid, dimname, dimlen, datamode) result(ncerr)
1967 
1968 !Arguments ------------------------------------
1969 !scalars
1970  integer,intent(in) :: ncid
1971  character(len=*),intent(in) :: dimname
1972  integer,intent(out) :: dimlen
1973  logical,optional,intent(in) :: datamode
1974 
1975 !Local variables-------------------------------
1976 !scalars
1977  integer :: dimid
1978 
1979 ! *********************************************************************
1980 
1981  ncerr = nf90_noerr
1982 
1983  if (present(datamode)) then
1984    if (datamode) then
1985      NCF_CHECK(nctk_set_datamode(ncid))
1986    end if
1987  end if
1988 
1989  ncerr = nf90_inq_dimid(ncid, dimname, dimid)
1990  if (ncerr == nf90_noerr) then
1991    ncerr = nf90_inquire_dimension(ncid, dimid, len=dimlen)
1992  end if
1993 
1994 end function nctk_get_dim

m_nctk/nctk_idgroup [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_idgroup

FUNCTION

  Return the nc identifier from the name of a group

SOURCE

346 integer function nctk_idgroup(ncid, grpname) result(grpid)
347 
348 !Arguments ------------------------------------
349  integer,intent(in) :: ncid
350  character(len=*),intent(in) :: grpname
351 
352 !Local variables-------------------------------
353 !scalars
354  integer :: ncerr
355  character(len=1000) :: msg
356 
357 ! *********************************************************************
358 
359 #ifdef HAVE_NETCDF
360  ncerr = nf90_inq_ncid(ncid, grpname, grpid)
361 
362  if (ncerr /= nf90_noerr) then
363    write(msg,'(6a)')&
364      "NetCDF library returned: `",trim(nf90_strerror(ncerr)), "`", ch10,&
365      "while trying to get the ncid of group: ",trim(grpname)
366    ABI_ERROR(msg)
367  end if
368 #else
369  ABI_ERROR("Netcdf support is not activated")
370  write(std_out,*)ncid,grpname
371 #endif
372 
373 end function nctk_idgroup

m_nctk/nctk_idname [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_idname

FUNCTION

  Return the nc identifier from the name of the variable

SOURCE

305 integer function nctk_idname(ncid, varname) result(varid)
306 
307 !Arguments ------------------------------------
308  integer,intent(in) :: ncid
309  character(len=*),intent(in) :: varname
310 
311 !Local variables-------------------------------
312 !scalars
313  integer :: ncerr
314  character(len=1000) :: msg
315 
316 ! *********************************************************************
317 
318 #ifdef HAVE_NETCDF
319  ncerr = nf90_inq_varid(ncid, varname, varid)
320 
321  if (ncerr /= nf90_noerr) then
322    write(msg,'(6a)')&
323      "NetCDF library returned: `",trim(nf90_strerror(ncerr)), "`", ch10,&
324      "while trying to get the ncid of variable: ",trim(varname)
325    ABI_ERROR(msg)
326  end if
327 #else
328  ABI_ERROR("Netcdf support is not activated")
329  write(std_out,*)ncid,varname
330 #endif
331 
332 end function nctk_idname

m_nctk/nctk_ncify [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_ncify

FUNCTION

  Append ".nc" to ipath if ipath does not end with ".nc"

SOURCE

386 function nctk_ncify(ipath) result(opath)
387 
388  character(len=*),intent(in) :: ipath
389  character(len=fnlen) :: opath
390 
391 ! *********************************************************************
392 
393  if (.not. endswith(ipath, ".nc")) then
394    opath = trim(ipath)//'.nc'
395  else
396    opath = ipath
397  end if
398 
399 end function nctk_ncify

m_nctk/nctk_open_create [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_open_create

FUNCTION

  Create and open the netcdf file. Return exis status.

INPUTS

  path=Name of the file
  comm=MPI communicator.

OUTPUT

  ncid=Netcdf identifier.

SOURCE

786 integer function nctk_open_create(ncid, path, comm) result(ncerr)
787 
788 !Arguments ------------------------------------
789  integer,intent(out) :: ncid
790  integer,intent(in) :: comm
791  character(len=*),intent(in) :: path
792 
793 !Local variables-------------------------------
794  integer :: input_len, cmode !, ii, ich
795  character(len=strlen) :: my_string
796 
797 ! *********************************************************************
798 
799  ! Always use mpiio mode (i.e. hdf5) if available so that one can perform parallel IO
800  if (nctk_has_mpiio) then
801    ncerr = nf90_einval
802 #ifdef HAVE_NETCDF_MPI
803    write(my_string,'(2a)') "- Creating HDf5 file with MPI-IO support: ",path
804    call wrtout(std_out,my_string)
805    ! Believe it or not, I have to use xmpi_comm_self even in sequential to avoid weird SIGSEV in the MPI layer!
806    ncerr = nf90_create(path, cmode=ior(ior(nf90_netcdf4, nf90_mpiio), nf90_write), ncid=ncid, &
807      comm=comm, info=xmpio_info)
808 #endif
809  else
810    ! Note that here we don't enforce nf90_netcdf4 hence the netcdf file with be in classic model.
811    write(my_string,'(2a)') "- Creating HDf5 file with MPI-IO support: ",path
812    call wrtout(std_out,my_string)
813    !ncerr = nf90_create(path, ior(nf90_clobber, nf90_write), ncid)
814    cmode = def_cmode_for_seq_create
815    ncerr = nf90_create(path, cmode=cmode, ncid=ncid)
816    if (xmpi_comm_size(comm) > 1) then
817      ABI_WARNING("netcdf without MPI-IO support with nprocs > 1!")
818    end if
819  end if
820  NCF_CHECK(ncerr)
821 
822  ! Write etsf_header: file format, version and conventions.
823  NCF_CHECK(nf90_put_att(ncid, NF90_GLOBAL, "file_format", etsfio_file_format))
824  NCF_CHECK(nf90_put_att(ncid, NF90_GLOBAL, "file_format_version", etsfio_version))
825  NCF_CHECK(nf90_put_att(ncid, NF90_GLOBAL, "Conventions", etsfio_conventions))
826 
827  ! Add info on the code that produced this file. These are extensions not in the standard.
828  NCF_CHECK(nf90_put_att(ncid, NF90_GLOBAL, "code", "Abinit"))
829  NCF_CHECK(nf90_put_att(ncid, NF90_GLOBAL, "abinit_version", abinit_version))
830 
831  ! Define the basic dimensions used in ETSF-IO files.
832  NCF_CHECK(nctk_def_basedims(ncid, defmode=.True.))
833 
834  if (len_trim(INPUT_STRING) /= 0) then
835    ! Write string with input.
836    my_string = trim(INPUT_STRING)
837    if (DTSET_IDX /= -1 .and. index(INPUT_STRING, "jdtset ") == 0) then
838      my_string = "jdtset " // trim(itoa(DTSET_IDX)) // "  " // trim(INPUT_STRING)
839    end if
840 
841    ! Since INPUT_STRING contains many control characters at the end (likely because it's a global var)
842    ! and we want to save space on disk, we cannot use trim_len and we have to find the last alphanum char in my_string.
843    input_len = len_trim(my_string)
844 #if 0
845    do ii=len(my_string), 1, -1
846      ich = iachar(my_string(ii:ii))
847      select case(ich)
848      case(0:32)  ! space, tab, or control character
849        !write(std_out, *)"space/tab/control at: ",ii, "iachar: ",iachar(my_string(ii:ii)), "char:", my_string(ii:ii)
850        cycle
851      case default
852        input_len = ii !; write(std_out, *)"Exiting at ii: ",ii, "with: ",my_string(ii:ii)
853        exit
854      end select
855    end do
856 #endif
857 
858    NCF_CHECK(nctk_def_dims(ncid, nctkdim_t("input_len", input_len)))
859    NCF_CHECK(nctk_def_arrays(ncid, nctkarr_t("input_string", "c", "input_len")))
860 
861    if (xmpi_comm_rank(comm) == 0) then
862      NCF_CHECK(nctk_set_datamode(ncid))
863      ! Pass my_string(1:input_len)) instead from trim(string) to avoid SIGSEV on higgs_intel_19.0_serial
864      NCF_CHECK(nf90_put_var(ncid, nctk_idname(ncid, "input_string"), my_string(1:input_len)))
865      NCF_CHECK(nctk_set_defmode(ncid))
866    end if
867  end if
868 
869 end function nctk_open_create

m_nctk/nctk_open_modify [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_open_modfy

FUNCTION

   Open an already existent netcdf file for modifications. Return exit status.

INPUTS

  path=File name.
  comm=MPI communicator.

OUTPUT

  ncid=Netcdf identifier.

SOURCE

890 integer function nctk_open_modify(ncid, path, comm) result(ncerr)
891 
892 !Arguments ------------------------------------
893  integer,intent(out) :: ncid
894  integer,intent(in) :: comm
895  character(len=*),intent(in) :: path
896 
897 ! *********************************************************************
898 
899  if (.not. nctk_has_mpiio .and. xmpi_comm_size(comm) > 1) then
900    ABI_ERROR("netcdf without MPI-IO support and nprocs > 1!")
901  end if
902 
903  if (xmpi_comm_size(comm) > 1 .or. nctk_has_mpiio) then
904    call wrtout(std_out, sjoin("- Opening HDf5 file with MPI-IO support:", path))
905 #ifdef HAVE_NETCDF_MPI
906    ncerr = nf90_open_par(path, cmode=ior(ior(nf90_netcdf4, nf90_mpiio), nf90_write), &
907      comm=comm, info=xmpio_info, ncid=ncid)
908    NCF_CHECK_MSG(ncerr, sjoin("nf90_open_par: ", path))
909 #else
910    ABI_ERROR("nprocs > 1 but netcdf does not support MPI-IO")
911 #endif
912  else
913    call wrtout(std_out, sjoin("- Opening netcdf file without MPI-IO support:", path))
914    ncerr = nf90_open(path, nf90_write, ncid)
915    NCF_CHECK_MSG(ncerr, sjoin("nf90_open: ", path))
916  end if
917 
918  ! Set file in define mode.
919  NCF_CHECK(nctk_set_defmode(ncid))
920 
921 end function nctk_open_modify

m_nctk/nctk_open_read [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_open_read

FUNCTION

   Open a netcdf file in read-only mode. Return exit status.

INPUTS

  ncid=Netcdf identifier.
  comm=MPI communicator.

SOURCE

729 integer function nctk_open_read(ncid, path, comm) result(ncerr)
730 
731 !Arguments ------------------------------------
732  integer,intent(out) :: ncid
733  integer,intent(in) :: comm
734  character(len=*),intent(in) :: path
735 
736 !Local variables-------------------------------
737  integer :: nprocs
738 
739 ! *********************************************************************
740  nprocs = xmpi_comm_size(comm)
741 
742  ! Enforce netcdf4 only if the communicator contains more than one processor.
743  if (nctk_has_mpiio .and. nprocs > 1) then
744 #ifdef HAVE_NETCDF_MPI
745    ncerr = nf90_open(path, mode=ior(ior(nf90_netcdf4, nf90_mpiio), nf90_nowrite),&
746                      comm=comm, info=xmpio_info, ncid=ncid)
747 #else
748    ncerr = nf90_einval
749    ABI_WARNING("Netcdf without MPI support. Cannot open file, will abort in caller")
750 #endif
751    NCF_CHECK_MSG(ncerr, sjoin("opening file:", path))
752  else
753    ncerr = nf90_open(path, mode=nf90_nowrite, ncid=ncid)
754    NCF_CHECK_MSG(ncerr, sjoin("opening file:", path))
755    !if (ncerr /= NC_EHDFERR) then
756    !  ncerr = nf90_open(path, mode=ior(ior(nf90_netcdf4), nf90_nowrite),&
757    !                    comm=comm, info=xmpio_info, ncid=ncid)
758    !end if
759    if (nprocs > 1) then
760      ncerr = nf90_einval
761      ABI_WARNING("netcdf without MPI-IO support with nprocs > 1! Will abort in the caller")
762    end if
763  endif
764 
765 end function nctk_open_read

m_nctk/nctk_read_datar [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

 nctk_read_datar

FUNCTION

  Read an array in real space in netcdf format

INPUTS

  path=Filename
  varname=Name of the variable to read.
  ngfft(18)=information about 3D FFT
  cplex=1 for real arrays (e.g. GS rhor), 2 for complex array.
  nfft=number of points in the real space FFT mesh treated by this MPI proc
  nspden=number of spin-density components
  comm_fft=MPI communicator (used only if MPI-FFT).
  fftn3_distrib(n3)=rank of the processors which own fft planes in 3rd dimension.
  ffti3_local(n3)=local index for 3d dimension
  datar(cplex*nfft,nspden)= array in real space.

OUTPUT

  Only writing

SOURCE

2226 integer function nctk_read_datar(path,varname,ngfft,cplex,nfft,nspden,&
2227    comm_fft,fftn3_distrib,ffti3_local,datar) result(ncerr)
2228 
2229 !Arguments ------------------------------------
2230 !scalars
2231  integer,intent(in) :: cplex,nfft,nspden,comm_fft
2232  character(len=*),intent(in) :: path,varname
2233 !arrays
2234  integer,intent(in) :: ngfft(18),fftn3_distrib(ngfft(3)),ffti3_local(ngfft(3))
2235  real(dp),intent(out) :: datar(cplex*nfft,nspden)
2236 
2237 !Local variables-------------------------------
2238 !scalars
2239  integer,parameter :: master=0
2240  integer :: ncid,varid,i3,nproc_fft,me_fft,i3_glob,n1,n2,n3,ispden
2241  logical :: ionode
2242 !arrays
2243  real(dp),allocatable :: glob_datar(:,:)
2244 
2245 ! *************************************************************************
2246 
2247  nproc_fft = xmpi_comm_size(comm_fft); me_fft = xmpi_comm_rank(comm_fft)
2248  n1 = ngfft(1); n2 = ngfft(2); n3 = ngfft(3)
2249 
2250  ! TODO: Be careful here because we should always create with HDF5 if available
2251  ! to avoid problems if we have to reread with nproc_fft > 1 and MPI-IO
2252  ionode = .True.
2253  if (nproc_fft == 1) then
2254    ncerr = nf90_open(path, mode=nf90_nowrite, ncid=ncid)
2255  else
2256    if (nctk_has_mpiio) then
2257      !write(std_out,*)"open_par: ",trim(path)
2258      !ncerr = nf90_open_par(path, nf90_nowrite,
2259      ! Don't know why but the format is not autodected!
2260      ncerr = nf90_einval
2261 #ifdef HAVE_NETCDF_MPI
2262      ncerr = nf90_open(path, mode=ior(ior(nf90_netcdf4, nf90_mpiio), nf90_nowrite),&
2263                        comm=comm_fft, info=xmpio_info, ncid=ncid)
2264 #endif
2265    else
2266      ! MPI-FFT without MPI-support. Only master does IO
2267      ionode = (me_fft == master); ncerr = nf90_noerr
2268      if (ionode) ncerr = nf90_open(path, nf90_nowrite, ncid)
2269    end if
2270  end if
2271  NCF_CHECK_MSG(ncerr, sjoin("opening file: ",path))
2272 
2273  NCF_CHECK(nf90_inq_varid(ncid, varname, varid))
2274  !write(std_out,*)"about to read varname, ngfft, cplex, nfft, nspden:", trim(varname), ngfft(:3), cplex,nfft,nspden
2275 
2276  ! netcdf array has shape [cplex, n1, n2, n3, nspden]
2277  if (nproc_fft == 1) then
2278    ! No MPI-FFT --> easy
2279    NCF_CHECK(nf90_get_var(ncid, varid, datar, start=[1,1,1,1], count=[cplex, n1, n2, n3, nspden]))
2280    NCF_CHECK(nf90_close(ncid))
2281 
2282  else
2283    ! Handle data distribution.
2284    ABI_CHECK(mod(ngfft(3), nproc_fft) == 0, "assuming mod(n3, nproc_fft) == 0")
2285 
2286    i3_glob = -1
2287    do i3=1,ngfft(3)
2288      if (fftn3_distrib(i3) == me_fft) then
2289         i3_glob = i3
2290         exit
2291      end if
2292    end do
2293    ABI_CHECK(i3_glob > 0, "negative i3_glob")
2294 
2295    if (nctk_has_mpiio) then
2296      ! Use parallel IO with collective calls.
2297      ncerr = nf90_einval
2298      NCF_CHECK(nctk_set_collective(ncid, varid))
2299 
2300      do ispden=1,nspden
2301        ncerr = nf90_get_var(ncid, varid, datar(:,ispden), start=[1,1,1,i3_glob,ispden], &
2302                 count=[cplex,ngfft(1),ngfft(2),ngfft(3)/nproc_fft,1])
2303        NCF_CHECK(ncerr)
2304      end do
2305    else
2306      ! MPI-FFT without MPI-IO. Master read and broadcast (requires more memory and communication)
2307      ABI_MALLOC(glob_datar, (cplex*product(ngfft(1:3)), nspden))
2308      if (ionode) then
2309        NCF_CHECK(nf90_get_var(ncid, varid, glob_datar, start=[1,1,1,1,1], count=[cplex,n1,n2,n3,nspden]))
2310      end if
2311 
2312      call distrib_datar(ngfft,cplex,nfft,nspden,glob_datar,master,comm_fft,fftn3_distrib,ffti3_local,datar)
2313      ABI_FREE(glob_datar)
2314    end if
2315 
2316    if (ionode) then
2317      NCF_CHECK(nf90_close(ncid))
2318    end if
2319  end if
2320 
2321 end function nctk_read_datar

m_nctk/nctk_set_atomic_units [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_set_atomic_units

FUNCTION

  Set the attributes "units" to "atomic units" and "scale_to_atomic_units" to one.

INPUTS

  ncid=Netcdf identifier.
  varname=Name of the variable

SOURCE

1256 integer function nctk_set_atomic_units(ncid, varname) result(ncerr)
1257 
1258 !Arguments ------------------------------------
1259  integer,intent(in) :: ncid
1260  character(len=*),intent(in) :: varname
1261 
1262 !Local variables-------------------------------
1263 !scalars
1264  integer :: varid
1265 
1266 ! *********************************************************************
1267 
1268  ncerr = nf90_noerr
1269 
1270  varid = nctk_idname(ncid, varname)
1271  NCF_CHECK(nf90_put_att(ncid, varid, "units", "atomic units"))
1272  NCF_CHECK(nf90_put_att(ncid, varid, "scale_to_atomic_units", one))
1273 
1274 end function nctk_set_atomic_units

m_nctk/nctk_set_collective [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_set_collective

FUNCTION

  Use collective IO for the given netcdf variable. Return exit status.

INPUTS

  ncid=Netcdf file identifier.
  varid=Netcdf variable identifier.

SOURCE

1099 integer function nctk_set_collective(ncid, varid) result(ncerr)
1100 
1101 !Arguments ------------------------------------
1102  integer,intent(in) :: ncid,varid
1103 
1104 ! *********************************************************************true
1105 
1106   ncerr = nf90_einval
1107 #ifdef HAVE_NETCDF_MPI
1108   ncerr = nf90_var_par_access(ncid, varid, nf90_collective)
1109 #else
1110   ABI_ERROR("nctk_set_collective should not be called if NETCDF does not support MPI-IO")
1111   ABI_UNUSED((/ncid, varid/))
1112 #endif
1113 
1114 end function nctk_set_collective

m_nctk/nctk_set_datamode [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_set_datamode

FUNCTION

  Set the file in data mode. Return exit status

INPUTS

  ncid=Netcdf identifier.
  [reserve]

OUTPUT

  ncerr=Exit status

SOURCE

1047 integer function nctk_set_datamode(ncid, reserve) result(ncerr)
1048 
1049 !Arguments ------------------------------------
1050  integer,intent(in) :: ncid
1051  logical,intent(in),optional :: reserve
1052 
1053 !Local variables-------------------------------
1054 !scalars
1055  logical :: do_reserve
1056 
1057 ! *********************************************************************
1058 
1059  do_reserve = .False.; if (present(reserve)) do_reserve = reserve
1060 
1061  ncerr = nf90_enddef(ncid)
1062 
1063  ! Use same trick as in etsf_io
1064  ! neded otherwise netcdf complains if the file is already in def mode.
1065  if (ncerr /= nf90_noerr .and. ncerr /= -38) then
1066    NCF_CHECK(ncerr)
1067  else
1068    ncerr = nf90_noerr
1069  end if
1070 
1071  return
1072 
1073  ! TODO
1074  if (do_reserve) then
1075    ncerr = nf90_enddef(ncid)
1076    !ncerr = nf90__enddef(ncid)
1077  else
1078    ncerr = nf90_enddef(ncid)
1079  end if
1080 
1081 end function nctk_set_datamode

m_nctk/nctk_set_default_for_seq [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_set_default_for_seq

FUNCTION

  Use netcdf classic mode for new files when only sequential-IO needs to be performed

SOURCE

281 subroutine nctk_use_classic_for_seq()
282 
283 ! *********************************************************************
284 
285 #ifdef HAVE_NETCDF
286  ! Use netcdf classic mode.
287  def_cmode_for_seq_create = ior(nf90_clobber, nf90_write)
288  ABI_COMMENT("Using netcdf-classic mode")
289 #endif
290 
291 end subroutine nctk_use_classic_for_seq

m_nctk/nctk_set_defmode [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_set_defmode

FUNCTION

   Set the file in define mode, return exit status.

INPUTS

  ncid=Netcdf identifier.

SOURCE

1011 integer function nctk_set_defmode(ncid) result(ncerr)
1012 
1013 !Arguments ------------------------------------
1014  integer,intent(in) :: ncid
1015 
1016 ! *********************************************************************
1017 
1018  ncerr = nf90_redef(ncid)
1019  ! Use same trick as in etsf_io
1020  if (ncerr /= nf90_noerr .and. ncerr /= -39) then
1021    NCF_CHECK(ncerr)
1022  else
1023    ncerr = nf90_noerr
1024  end if
1025 
1026 end function nctk_set_defmode

m_nctk/nctk_string_from_occopt [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_string_from_occopt

FUNCTION

SOURCE

413 pure function nctk_string_from_occopt(occopt) result(smearing)
414 
415  integer,intent(in) :: occopt
416  character(len=etsfio_charlen) :: smearing
417 
418 ! *********************************************************************
419 
420  select case (occopt)
421  case (3)
422    smearing = "Fermi-Dirac"
423  case (4)
424    smearing = "cold smearing of N. Marzari with minimization of the bump"
425  case (5)
426    smearing = "cold smearing of N. Marzari with monotonic function in the tail"
427  case (6)
428    smearing = "Methfessel and Paxton"
429  case (7)
430    smearing = "gaussian"
431  case (8)
432    smearing = "uniform"
433  case default
434    smearing = "none"
435  end select
436 
437 end function nctk_string_from_occopt

m_nctk/nctk_test_mpiio [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_test_mpiio

FUNCTION

  Test at run-time whether the netcdf library supports parallel IO and
  set the value of the module variable `nctk_has_mpiio`.
  This is a COLLECTIVE routine that should be called by all processors
  in MPI_COMM_WORLD at the beginning of the calculation

INPUTS

  [print_warning]=TRUE if a warning about paral_kgb use has to be printed
                  Optional, default=yes

SOURCE

565 subroutine nctk_test_mpiio(print_warning)
566 
567  logical,intent(in),optional :: print_warning
568 
569 !Local variables-------------------------------
570 !scalars
571  logical :: my_print_warning
572 #ifdef HAVE_NETCDF_MPI
573  integer,parameter :: master=0
574  integer :: ierr,ncid,ncerr
575  character(len=500) :: msg
576  character(len=fnlen) :: apath
577 #endif
578 
579 ! *********************************************************************
580 
581  nctk_has_mpiio = .False.
582  my_print_warning=.true. ; if (present(print_warning)) my_print_warning=print_warning
583 
584 !FIXME nf90create fails when using NVHPC
585 ! This might be due to my environment, maybe not, need to investigate this...
586 !!#ifndef FC_NVHPC
587 #ifdef HAVE_NETCDF_MPI
588  if (xmpi_comm_rank(xmpi_world) == master) then
589    ! Try to open a file with hdf5.
590    apath = pick_aname()
591    ncerr = nf90_create(apath, cmode=ior(ior(nf90_netcdf4, nf90_mpiio), nf90_write), ncid=ncid, &
592      comm=xmpi_comm_self, info=xmpio_info)
593 
594    if (ncerr == nf90_noerr) then
595      nctk_has_mpiio = .True.
596      call wrtout(std_out," Netcdf library supports MPI-IO", "COLL")
597    else if (ncerr == nf90_enopar) then
598      ! This is the value returned by the C function ifndef USE_PARALLEL
599      ABI_WARNING(sjoin("Netcdf lib does not support MPI-IO and: ", nf90_strerror(ncerr)))
600      nctk_has_mpiio = .False.
601    else
602      ! Maybe something wrong in the low-level layer!
603      ABI_WARNING(sjoin("Strange, netcdf seems to support MPI-IO but: ", nf90_strerror(ncerr)))
604      nctk_has_mpiio = .False.
605    end if
606 
607    ncerr = nf90_close(ncid)
608    call delete_file(apath, ierr)
609  end if
610 
611  ! Master broadcast nctk_has_mpiio
612  call xmpi_bcast(nctk_has_mpiio,master,xmpi_world,ierr)
613 
614  if ((.not.nctk_has_mpiio).and.my_print_warning) then
615    write(msg,"(5a)") &
616       "The netcdf library does not support parallel IO, see message above",ch10,&
617       "Abinit won't be able to produce files in parallel e.g. when paral_kgb==1 is used.",ch10,&
618       "Action: install a netcdf4+HDF5 library with MPI-IO support."
619    ABI_WARNING(msg)
620  end if
621 #endif
622 !!#endif
623 
624 #ifdef HAVE_NETCDF_DEFAULT
625  if (.not. nctk_has_mpiio) then
626    ABI_ERROR("--netcdf-default is on but netcdf library does not support MPI-IO. Aborting now")
627  end if
628 #endif
629 
630 end subroutine nctk_test_mpiio

m_nctk/nctk_try_fort_or_ncfile [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_try_fort_or_ncfile

FUNCTION

  If filename does not exist, a similar file with extension `.nc` is tried
  This trick is used to run the Abinit test suite in netcdf mode without changing the input files.
  The modification (if any) is logged to unit (Default: std_out)

SIDE EFFECTS

  filename=Tentative filename in input. Changed to netcdf file if input filename does not exist
   and a file with netcdf extension is found.

OUTPUT

  errmsg=String with error message if return value /= 0

SOURCE

516 integer function nctk_try_fort_or_ncfile(filename, errmsg, unit) result(ierr)
517 
518 !Arguments ------------------------------------
519  character(len=*),intent(inout) :: filename
520  character(len=*),intent(out) :: errmsg
521  integer,optional,intent(in) :: unit
522 
523 !Local variables-------------------------------
524  integer :: unt
525 
526 ! *********************************************************************
527 
528  unt = std_out; if (present(unit)) unt = unit
529  ierr = 0; errmsg = ""
530 
531  if (.not.file_exists(filename)) then
532    ! Try netcdf exists.
533    if (file_exists(nctk_ncify(filename))) then
534      if (unt /= dev_null) then
535        write(unt,"(3a)")"- File: ",trim(filename)," does not exist but found netcdf file with similar name."
536      end if
537      filename = nctk_ncify(filename)
538    end if
539    if (.not. file_exists(filename)) then
540      ierr = 1; errmsg = 'Cannot find file: '//trim(filename)
541    end if
542  end if
543 
544 end function nctk_try_fort_or_ncfile

m_nctk/nctk_write_datar [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

 nctk_write_datar

FUNCTION

  Write an array in real space in netcdf format

INPUTS

  path=Filename
  varname=Name of the variable to write.
  ngfft(18)=information about 3D FFT
  cplex=1 for real arrays (e.g. GS rhor), 2 for complex array.
  nfft=number of points in the real space FFT mesh treated by this MPI proc
  nspden=number of spin-density components
  comm_fft=MPI communicator (used only if MPI-FFT).
  fftn3_distrib(n3)=rank of the processors which own fft planes in 3rd dimension.
  ffti3_local(n3)=local index for 3d dimension
  datar(cplex*nfft,nspden)= array in real space.
  [action]

OUTPUT

  Only writing

SOURCE

2024 integer function nctk_write_datar(varname,path,ngfft,cplex,nfft,nspden,&
2025    comm_fft,fftn3_distrib,ffti3_local,datar,action) result(ncerr)
2026 
2027 !Arguments ------------------------------------
2028 !scalars
2029  integer,intent(in) :: cplex,nfft,nspden,comm_fft
2030  character(len=*),intent(in) :: path,varname
2031  character(len=*),optional,intent(in) :: action
2032 !arrays
2033  integer,intent(in) :: ngfft(18),fftn3_distrib(ngfft(3)),ffti3_local(ngfft(3))
2034  real(dp),target,intent(in) :: datar(cplex*nfft,nspden)
2035 
2036 !Local variables-------------------------------
2037 !scalars
2038  integer,parameter :: master=0
2039  integer :: ncid,varid,i3,nproc_fft,me_fft,i3_glob,n1,n2,n3,ispden,cmode
2040  logical :: ionode
2041  character(len=nctk_slen) :: cplex_name,my_action
2042  !character(len=500) :: msg
2043 !arrays
2044  real(dp),allocatable :: glob_datar(:,:)
2045 
2046 ! *************************************************************************
2047 
2048  ! FIXME: Default should be open but this enters into conflict with the abi_estf stuff!
2049  ! if we are using MPI-IO since file is not open with HDF5.
2050  !my_action = "create"; if (present(action)) my_action = action
2051  my_action = "open"; if (present(action)) my_action = action
2052 
2053  nproc_fft = xmpi_comm_size(comm_fft); me_fft = xmpi_comm_rank(comm_fft)
2054  n1 = ngfft(1); n2 = ngfft(2); n3 = ngfft(3)
2055 
2056  ! TODO: Be careful here because we should always create with HDF5 if available
2057  ! to avoid problems if we have to reread with nproc_fft > 1 and MPI-IO
2058  ionode = .True.; ncerr = nf90_noerr
2059  if (nproc_fft == 1) then
2060    select case(my_action)
2061    case ("open")
2062      ncerr = nf90_open(path, mode=nf90_write, ncid=ncid)
2063    case ("create")
2064      ncerr = nctk_open_create(ncid, path, comm_fft)
2065    case default
2066      ABI_ERROR(sjoin("Wrong action: ", my_action))
2067    end select
2068 
2069  else
2070    if (nctk_has_mpiio) then
2071      call wrtout(std_out, strcat("nctk_write_datar: using MPI-IO to write ", varname, path), "COLL")
2072 
2073      ncerr = nf90_einval
2074 #ifdef HAVE_NETCDF_MPI
2075      select case(my_action)
2076      case ("open")
2077        ncerr = nf90_open(path, mode=nf90_write, comm=comm_fft, info=xmpio_info, ncid=ncid)
2078      case ("create")
2079        ncerr = nf90_create(path, cmode=ior(ior(nf90_netcdf4, nf90_mpiio), nf90_write), &
2080          comm=comm_fft, info=xmpio_info, ncid=ncid)
2081      case default
2082        ABI_ERROR(strcat("Wrong action:", my_action))
2083      end select
2084 #endif
2085    else
2086      ! MPI-FFT without MPI-support. Only master performs IO
2087      ionode = (me_fft == master)
2088      if (ionode) then
2089        select case(my_action)
2090        case ("open")
2091          ncerr = nf90_open(path, mode=nf90_write, ncid=ncid)
2092        case ("create")
2093          !ncerr = nf90_create(path, cmode=nf90_clobber, ncid=ncid)
2094          cmode = def_cmode_for_seq_create
2095          ncerr = nf90_create(path, cmode=cmode, ncid=ncid)
2096        case default
2097          ABI_ERROR(strcat("Wrong action:", my_action))
2098        end select
2099      end if
2100    end if
2101  end if
2102  NCF_CHECK_MSG(ncerr, sjoin("opening file:", path))
2103 
2104  if (ionode) then
2105    ! Define dims and variables
2106    !write(std_out,*)"defing dims",trim(varname)," in file: ",path
2107    cplex_name = strcat("real_or_complex_", varname)
2108    ncerr = nctk_def_dims(ncid, [&
2109      nctkdim_t(cplex_name, cplex),&
2110      nctkdim_t("number_of_grid_points_vector1", n1),&
2111      nctkdim_t("number_of_grid_points_vector2", n2),&
2112      nctkdim_t("number_of_grid_points_vector3", n3),&
2113      nctkdim_t("number_of_components", nspden)], defmode=.True.)
2114    NCF_CHECK(ncerr)
2115 
2116    ncerr = nctk_def_one_array(ncid, nctkarr_t(name=varname, dtype="dp", shape_str=strcat(cplex_name, &
2117 ", number_of_grid_points_vector1, number_of_grid_points_vector2, number_of_grid_points_vector3, number_of_components")),&
2118    varid=varid)
2119 
2120    ! Add attributes
2121    varid = nctk_idname(ncid, varname)
2122    NCF_CHECK(nf90_put_att(ncid, varid, "units", "atomic units"))
2123    NCF_CHECK(nf90_put_att(ncid, varid, "scale_to_atomic_units", one))
2124  end if
2125 
2126  if (nproc_fft == 1) then
2127    ! no MPI-FFT --> write data directly.
2128    varid = nctk_idname(ncid, varname)
2129    NCF_CHECK(nctk_set_datamode(ncid))
2130    NCF_CHECK(nf90_put_var(ncid, varid, datar, start=[1,1,1,1,1], count=[cplex, n1, n2, n3, nspden]))
2131    NCF_CHECK(nf90_close(ncid))
2132 
2133  else
2134    ! Must handle data distribution.
2135    ABI_CHECK(mod(n3, nproc_fft) == 0, "assuming mod(n3, nproc_fft) == 0")
2136 
2137    i3_glob = -1
2138    do i3=1,ngfft(3)
2139      if (fftn3_distrib(i3) == me_fft) then
2140         i3_glob = i3
2141         exit
2142      end if
2143    end do
2144    ABI_CHECK(i3_glob > 0, "negative i3_glob")
2145 
2146    !do i3=i3_glob,ngfft(3)
2147    !  if (fftn3_distrib(i3) /= me_fft) exit
2148    !  !assert all(ffn3_distrib(i3_glob:i3_glob -1 + ngfft(3) / nproc_fft) == me_fft)
2149    !end do
2150    !i3_glob = i3- 1
2151    !print*,"i3_glob",i3_glob
2152 
2153    if (ionode) then
2154      NCF_CHECK(nf90_enddef(ncid))
2155    end if
2156 
2157    ! Array on disk has shape [cplex, n1, n2, n3, nspden]
2158    if (nctk_has_mpiio) then
2159      ! Use collective IO.
2160      ncerr = nf90_einval
2161      NCF_CHECK(nctk_set_collective(ncid, varid))
2162 
2163      do ispden=1,nspden
2164        ncerr = nf90_put_var(ncid, varid, datar(:,ispden), start=[1,1,1,i3_glob,ispden], &
2165                 count=[cplex,n1,n2,n3/nproc_fft,1])
2166        NCF_CHECK(ncerr)
2167      end do
2168    else
2169      ! MPI-FFT without MPI-IO. Collect data (requires more memory and communication)
2170      ABI_MALLOC(glob_datar, (cplex*product(ngfft(1:3)), nspden))
2171      call collect_datar(ngfft,cplex,nfft,nspden,datar,comm_fft,fftn3_distrib,ffti3_local,glob_datar,master=master)
2172 
2173      if (ionode) then
2174        ! Write global array.
2175        NCF_CHECK(nf90_put_var(ncid, varid, glob_datar, start=[1,1,1,1,1], count=[cplex,n1,n2,n3,nspden]))
2176      end if
2177      ABI_FREE(glob_datar)
2178    end if
2179 
2180    if (ionode) then
2181      NCF_CHECK(nf90_close(ncid))
2182    end if
2183  end if
2184 
2185  !ok = .True.
2186  ! Sequential IO
2187  !do rank=0,nproc_fft-1
2188  !  if (rank == me_fft) then
2189  !     ncerr = nf90_open(path, mode=nf90_write, ncid=ncid)
2190  !     do ispden=1,nspden
2191  !       ncerr = nf90_put_var(ncid, varid, datar(1:,ispden), start=[1,1,1,i3_glob,ispden], &
2192  !                count=[cplex,ngfft(1),ngfft(2),ngfft(3)/nproc_fft,1])
2193  !     end do
2194  !     ncerr = nf90_close(ncid)
2195  !  end if
2196  !  call xmpi_barrier(comm_fft)
2197  !end do
2198 
2199 end function nctk_write_datar

m_nctk/nctk_write_dpscalars [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_write_dpscalars

FUNCTION

  Write a list of **scalar** real(dp) variables. Return immediately if error

INPUTS

  ncid=Netcdf identifier.
  varnames(:)=List of strings with the name of the variables
  values(:)=List of real(dp) scalars
  [datamode]=If True, the nc file is set in data mode (default=False)

SOURCE

1766 integer function nctk_write_dpscalars(ncid, varnames, values, datamode) result(ncerr)
1767 
1768 !Arguments ------------------------------------
1769 !scalars
1770  integer,intent(in) :: ncid
1771  logical,optional,intent(in) :: datamode
1772 !arrays
1773  real(dp),intent(in) :: values(:)
1774  character(len=*),intent(in) :: varnames(:)
1775 
1776 !Local variables-------------------------------
1777 !scalars
1778  integer :: ii,varid
1779 
1780 ! *********************************************************************
1781 
1782  ncerr = nf90_noerr
1783 
1784  ABI_CHECK(size(varnames) == size(values), "Different size in varnames, values")
1785 
1786  if (present(datamode)) then
1787    if (datamode) then
1788      NCF_CHECK(nctk_set_datamode(ncid))
1789    end if
1790  end if
1791 
1792  do ii=1,size(varnames)
1793    NCF_CHECK(nf90_inq_varid(ncid, varnames(ii), varid))
1794    NCF_CHECK(nf90_put_var(ncid, varid, values(ii)))
1795  end do
1796 
1797 end function nctk_write_dpscalars

m_nctk/nctk_write_ibz [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_write_ibz

FUNCTION

  Write the list of the k-points in the IBZ with the corresponding weights in Netcdf format.
  Mainly used for passing data to AbiPy. This routine should be called by master only.

INPUTS

  fname=File name
  kpoints(:,:)=List of k-points
  weights(:)=K-point weights

OUTPUT

  ncerr=Exit status

SOURCE

1909 integer function nctk_write_ibz(fname, kpoints, weights) result(ncerr)
1910 
1911 !Arguments ------------------------------------
1912 !scalars
1913  character(len=*),intent(in) :: fname
1914 !arrays
1915  real(dp),intent(in) :: kpoints(:,:),weights(:)
1916 
1917 !Local variables-------------------------------
1918 !scalars
1919  integer :: nkpts,ncid
1920 
1921 ! *********************************************************************
1922 
1923  ABI_CHECK(size(kpoints, dim=2) == size(weights), "size(kpoints, dim=2) != size(weights)")
1924  nkpts = size(kpoints, dim=2)
1925 
1926  NCF_CHECK_MSG(nctk_open_create(ncid, fname, xmpi_comm_self), sjoin("Creating:", fname))
1927 
1928  ncerr = nctk_def_dims(ncid, [ &
1929    nctkdim_t("number_of_reduced_dimensions",3), nctkdim_t("number_of_kpoints", nkpts)], defmode=.True.)
1930  NCF_CHECK(ncerr)
1931 
1932  ncerr = nctk_def_array_list(ncid, [&
1933    nctkarr_t('reduced_coordinates_of_kpoints', "dp", "number_of_reduced_dimensions, number_of_kpoints"),&
1934    nctkarr_t('kpoint_weights', "dp", "number_of_kpoints")])
1935  NCF_CHECK(ncerr)
1936 
1937  NCF_CHECK(nctk_set_datamode(ncid))
1938 
1939  ncerr = nf90_put_var(ncid, nctk_idname(ncid, 'reduced_coordinates_of_kpoints'), kpoints)
1940  NCF_CHECK(ncerr)
1941 
1942  NCF_CHECK(nf90_put_var(ncid, nctk_idname(ncid, 'kpoint_weights'), weights))
1943 
1944  NCF_CHECK(nf90_close(ncid))
1945 
1946 end function nctk_write_ibz

m_nctk/nctk_write_iscalars [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  nctk_write_iscalars

FUNCTION

  Write a list of **scalar** integer variables. Return immediately if error

INPUTS

  ncid=Netcdf identifier.
  varnames(:)=List of strings with the name of the variables
  values(:)=List of integer scalars
  [datamode]=If True, the nc file is set in data mode (default=False)

OUTPUT

  ncerr=Exit status

SOURCE

1718 integer function nctk_write_iscalars(ncid, varnames, values, datamode) result(ncerr)
1719 
1720 !Arguments ------------------------------------
1721 !scalars
1722  integer,intent(in) :: ncid
1723  logical,optional,intent(in) :: datamode
1724 !arrays
1725  integer,intent(in) :: values(:)
1726  character(len=*),intent(in) :: varnames(:)
1727 
1728 !Local variables-------------------------------
1729 !scalars
1730  integer :: ii,varid
1731 
1732 ! *********************************************************************
1733 
1734  ABI_CHECK(size(varnames) == size(values), "Different size in varnames, values")
1735 
1736  ncerr = nf90_noerr
1737  if (present(datamode)) then
1738    if (datamode) then
1739      NCF_CHECK(nctk_set_datamode(ncid))
1740    end if
1741  end if
1742 
1743  do ii=1,size(varnames)
1744    NCF_CHECK_MSG(nf90_inq_varid(ncid, varnames(ii), varid), sjoin("Inquiring: ", varnames(ii)))
1745    NCF_CHECK(nf90_put_var(ncid, varid, values(ii)))
1746  end do
1747 
1748 end function nctk_write_iscalars

m_nctk/nctkarr_t [ Types ]

[ Top ] [ m_nctk ] [ Types ]

NAME

 nctkarr_t

FUNCTION

  Stores the name and the shape of a netcdf array

SOURCE

135  type,public :: nctkarr_t
136    character(len=nctk_slen) :: name        ! name of the array.
137    character(len=4) :: dtype               ! string specifying the type.
138    character(len=nctk_slen) :: shape_str   ! string with the shape. e.g. "dim1, dim2" for [dim1, dim2] array.
139  end type nctkarr_t

m_nctk/nctkerr_t [ Types ]

[ Top ] [ m_nctk ] [ Types ]

NAME

 nctkerr_t

FUNCTION

SOURCE

 96  type,private :: nctkerr_t
 97 #ifdef HAVE_NETCDF
 98    integer :: ncerr = nf90_noerr
 99 #else
100    integer :: ncerr = 0
101 #endif
102    integer :: line = 0
103    character(len=fnlen) :: file = "Dummy File"
104    character(len=2048) :: msg="No error detected"
105  end type nctkerr_t

m_nctk/nctkvar_t [ Structures ]

[ Top ] [ m_nctk ] [ Structures ]

NAME

  nctkvar_t

FUNCTION

  This structure stores variable information, such as
  name, NetCDF id, type, shape and dimensions. It contains the following elements:

SOURCE

152  type nctkvar_t
153 
154    integer :: id
155    ! the id used by NetCDF to access this variable.
156 
157    integer :: xtype
158    ! the type of the variable
159 
160    integer :: ndims
161    ! the number of dimensions (0 for scalar variable).
162 
163    integer :: natts
164    ! The number of attributes associated to the variable
165 
166    character(len=nctk_slen) :: name
167    ! the variable name.
168 
169    character(len=nctk_slen) :: dimnames(nctk_max_dims)
170    ! the name corresponding to each dimension
171    ! Only if array variable, only (1:ndims) are relevent
172 
173    integer :: dimids(nctk_max_dims) = -1
174    ! The id of the dimensions. only (1:ndims) are relevent
175 
176    integer :: dimlens(nctk_max_dims) = 0
177    ! the size for each dimension if array variable, only (1:ndims) are relevent
178 
179    !character(len=nctk_slen) :: dimnames(nctk_max_dims)
180    !character(len=nctk_slen), pointer :: ncattnames(:)
181    ! * ncattnames: the name corresponding to all associated attributes
182 
183  end type nctkvar_t

m_nctk/str2xtype [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  str2xtype

FUNCTION

  Return the netcdf type from a string. Possible values:
    c or ch   for NF90_CHAR
    i or int  for NF90_INTtrue
   sp         for NF90_FLOAT
   dp         for NF90_DOUBLE

SOURCE

648 integer function str2xtype(string) result(xtype)
649 
650 !Arguments ------------------------------------
651  character(len=*),intent(in) :: string
652 
653 ! *********************************************************************
654 
655  !Type  FORTRAN API Mnemonic    Bits
656  !byte      NF90_BYTE           8
657  !char      NF90_CHAR           8
658  !short     NF90_SHORT          16
659  !int       NF90_INT            32
660  !float     NF90_FLOAT          32
661  !double    NF90_DOUBLE         64
662 
663  select case (string)
664  case ("c", "ch", "char")
665    xtype = NF90_CHAR
666  case ("i", "int")
667    xtype = NF90_INT
668  case ("sp")
669    xtype = NF90_FLOAT
670  case ("dp")
671    xtype = NF90_DOUBLE
672  case default
673    ABI_ERROR(sjoin("Invalid string type:", string))
674  end select
675 
676 end function str2xtype

m_nctk/var_from_id [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  var_from_id

FUNCTION

  Initialize a nctkvar_t object from the variable id

INPUTS

  ncid=NC file handle
  varid=Variable ID

OUTPUT

  var<nctkvar_t>=Info on the variable.

SOURCE

2488 subroutine var_from_id(ncid, varid, var)
2489 
2490 !Arguments ------------------------------------
2491  integer, intent(in) :: ncid, varid
2492  type(nctkvar_t), intent(out) :: var
2493 
2494 !Local variables-------------------------------
2495 !scalars
2496  integer :: ii, ncerr
2497  !character(len=NF90_MAX_NAME) :: ncname
2498 
2499 ! *********************************************************************
2500 
2501  ! Get info about the variable.
2502  var%id = varid
2503  ncerr = nf90_inquire_variable(ncid, var%id, &
2504    name=var%name, xtype=var%xtype, ndims=var%ndims, dimids=var%dimids, natts=var%natts)
2505  NCF_CHECK(ncerr)
2506 
2507  ! Get info about dimensions.
2508  if (var%ndims > 0) then
2509    do ii=1,var%ndims
2510      ncerr = nf90_inquire_dimension(ncid, var%dimids(ii), len=var%dimlens(ii), name=var%dimnames(ii))
2511      NCF_CHECK(ncerr)
2512    end do
2513  end if
2514 
2515  ! Get the number of attributes and their names.
2516  !if (var%natts > 0) then
2517  !   do ii=1,var%natts
2518  !      ncerr = nf90_inq_attname(ncid, var%id, ii, var%attnames(ii))
2519  !   end do
2520  !end if
2521 
2522 end subroutine var_from_id

m_nctk/var_from_name [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

  var_from_name

FUNCTION

  Initialize a nctkvar_t object from the variable name

INPUTS

  ncid=NC file handle
  name=Variable name

OUTPUT

  var<nctkvar_t>=Info on the variable.

SOURCE

2543 subroutine var_from_name(ncid, name, var)
2544 
2545 !Arguments ------------------------------------
2546  integer, intent(in) :: ncid
2547  character(len=*),intent(in) :: name
2548  type(nctkvar_t), intent(out) :: var
2549 
2550 !Local variables-------------------------------
2551 !scalars
2552  integer :: varid
2553 
2554 ! *********************************************************************
2555 
2556  varid = nctk_idname(ncid, name)
2557  call var_from_id(ncid, varid, var)
2558 
2559 end subroutine var_from_name

m_nctk/write_var_netcdf [ Functions ]

[ Top ] [ m_nctk ] [ Functions ]

NAME

 write_var_netcdf

FUNCTION

 Write variable into a netcdf dataset

 TODO:
  Remove!

INPUTS

 arr_int
 arr_real
 marr
 narr
 typevar
 varname

OUTPUT

  (only writing)

SOURCE

2815 subroutine write_var_netcdf(arr_int,arr_real,marr,narr,ncid,typevar,varname)
2816 
2817 !Arguments ------------------------------------
2818 !scalars
2819  integer,intent(in) :: narr,marr,ncid
2820  character(len=*),intent(in) :: varname
2821  character(len=3),intent(in) :: typevar
2822 !arrays
2823  integer,intent(in) :: arr_int(marr)
2824  real(dp),intent(in) :: arr_real(marr)
2825 
2826 !Local variables-------------------------------
2827 !scalars
2828  integer :: var_id,var_type,vardim_id,ncerr
2829  !character(len=500) :: msg
2830 
2831 ! *************************************************************************
2832 
2833  !write(std_out,*)"about to write varname: ",trim(varname)
2834 
2835 #if defined HAVE_NETCDF
2836  if (ncid>0) then
2837 !  ### Put the file in definition mode
2838    ncerr=nf90_redef(ncid)
2839    if (ncerr/=NF90_NOERR.and.ncerr/=NF90_EINDEFINE) then
2840      NCF_CHECK_MSG(ncerr,'nf90_redef')
2841    end if
2842 !  ### Define the dimensions
2843    if (narr==1)then
2844      ncerr=nf90_inq_dimid(ncid,'one',vardim_id)
2845      NCF_CHECK_MSG(ncerr,'nf90_inq_varid')
2846    else
2847      ncerr=nf90_def_dim(ncid,trim(varname),narr,vardim_id)
2848      NCF_CHECK_MSG(ncerr,'nf90_def_dim')
2849    end if
2850 !  ### Define the variables
2851    if (typevar=='INT') then
2852      var_type=NF90_INT
2853    else if (typevar=='DPR') then
2854      var_type=NF90_DOUBLE
2855    end if
2856    ncerr=nf90_def_var(ncid, trim(varname), var_type, vardim_id, var_id)
2857    NCF_CHECK_MSG(ncerr,'nf90_def_var')
2858 !  ### Put the file in data mode
2859    ncerr=nf90_enddef(ncid)
2860    if (ncerr/=NF90_NOERR.and.ncerr/=NF90_ENOTINDEFINE) then
2861      NCF_CHECK_MSG(ncerr,'nf90_enddef')
2862    end if
2863 !  ### Write variables into the dataset
2864    if (typevar=='INT') then
2865      ncerr=nf90_put_var(ncid,var_id,arr_int,start=(/1/),count=(/narr/))
2866    else if (typevar=='DPR') then
2867      ncerr=nf90_put_var(ncid,var_id,arr_real,start=(/1/),count=(/narr/))
2868    end if
2869    NCF_CHECK_MSG(ncerr,'nf90_put_var')
2870  end if
2871 #endif
2872 
2873 end subroutine write_var_netcdf