TABLE OF CONTENTS


ABINIT/m_dfpt_io [ Modules ]

[ Top ] [ Modules ]

NAME

 m_dfpt_io

FUNCTION

 Module containing the methods used to do IO on DFPT results. 

COPYRIGHT

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

PARENTS

CHILDREN

SOURCE

21 #if defined HAVE_CONFIG_H
22 #include "config.h"
23 #endif
24 
25 #include "abi_common.h"
26 
27 MODULE m_dfpt_io
28 
29  use defs_basis
30  use m_profiling_abi
31  use m_errors
32  use m_nctk
33 #ifdef HAVE_NETCDF
34  use netcdf
35 #endif
36 
37  implicit none
38 
39  private 

m_dfpt_io/elast_ncwrite [ Functions ]

[ Top ] [ m_dfpt_io ] [ Functions ]

NAME

  elast_ncwrite

FUNCTION

 Writes the elastic constants tensors (clamped-ion, relaxed-ion with and
 without stress corrections) to a netCDF file.

INPUTS

 compl=relaxed-ion compliance tensor(without stress correction) (6*6) in
       Voigt notation
 compl_clamped=clamped-ion compliance tensor(without stress correction) (6*6)
               in Voigt notation
 compl_stress=relaxed-ion compliance tensor(with stress correction) (6*6) in
              Voigt notation
 elast=relaxed-ion elastic tensor(without stress correction) (6*6) in
       Voigt notation
 elast_clamped=clamped-ion elastic tensor(without stress correction) (6*6)
               in Voigt notation
 elast_stress=relaxed-ion elastic tensor(with stress correction) (6*6) in
              Voigt notation
 ncid=NC file handle (open in the caller)

OUTPUT

 Only writing

NOTES

 Units are GPa for elastic constants and GPa^-1 for compliance constants

PARENTS

      anaddb

CHILDREN

SOURCE

 83 subroutine elast_ncwrite(compl,compl_clamped,compl_stress,elast,elast_clamped,elast_stress,ncid)
 84 
 85 
 86 !This section has been created automatically by the script Abilint (TD).
 87 !Do not modify the following lines by hand.
 88 #undef ABI_FUNC
 89 #define ABI_FUNC 'elast_ncwrite'
 90 !End of the abilint section
 91 
 92  implicit none
 93 
 94 !Arguments ------------------------------------
 95 !scalars
 96  integer,intent(in) :: ncid 
 97 !arrays
 98  real(dp),intent(in) :: compl(6,6), compl_clamped(6,6),compl_stress(6,6)
 99  real(dp),intent(in) :: elast(6,6), elast_clamped(6,6),elast_stress(6,6)
100 
101 !Local variables-------------------------------
102 !scalars
103 #ifdef HAVE_NETCDF
104  integer :: ncerr
105 
106 ! *************************************************************************
107 
108 ! Define dimensions
109  NCF_CHECK(nctk_def_basedims(ncid, defmode=.True.))
110 
111 !arrays
112  ncerr = nctk_def_arrays(ncid, [&
113    nctkarr_t('compliance_constants_relaxed_ion', "dp", 'six, six'), &
114    nctkarr_t('compliance_constants_clamped_ion', "dp", 'six, six'), &
115    nctkarr_t('compliance_constants_relaxed_ion_stress_corrected', "dp", "six, six"), &
116    nctkarr_t('elastic_constants_relaxed_ion', "dp", 'six, six'), &
117    nctkarr_t('elastic_constants_clamped_ion', "dp", 'six, six'), &
118    nctkarr_t('elastic_constants_relaxed_ion_stress_corrected', "dp", 'six, six')])
119  NCF_CHECK(ncerr)
120 
121  ! Write variables. 
122  NCF_CHECK(nctk_set_datamode(ncid))
123  NCF_CHECK(nf90_put_var(ncid, vid('compliance_constants_relaxed_ion'), compl))
124  NCF_CHECK(nf90_put_var(ncid, vid('compliance_constants_clamped_ion'), compl_clamped))
125  NCF_CHECK(nf90_put_var(ncid, vid('compliance_constants_relaxed_ion_stress_corrected'), compl_stress))
126  NCF_CHECK(nf90_put_var(ncid, vid('elastic_constants_relaxed_ion'), elast))
127  NCF_CHECK(nf90_put_var(ncid, vid('elastic_constants_clamped_ion'), elast_clamped))
128  NCF_CHECK(nf90_put_var(ncid, vid('elastic_constants_relaxed_ion_stress_corrected'), elast_stress))
129 
130 #else
131  MSG_ERROR("Netcdf support not enabled")
132  ABI_UNUSED((/ncid/))
133 #endif
134 
135 contains
136  integer function vid(vname) 
137 
138 
139 !This section has been created automatically by the script Abilint (TD).
140 !Do not modify the following lines by hand.
141 #undef ABI_FUNC
142 #define ABI_FUNC 'vid'
143 !End of the abilint section
144 
145    character(len=*),intent(in) :: vname
146    vid = nctk_idname(ncid, vname)
147  end function vid
148 
149 end subroutine elast_ncwrite