TABLE OF CONTENTS


ABINIT/outddbnc [ Functions ]

[ Top ] [ Functions ]

NAME

 outddbnc

FUNCTION

 Write Derivative DataBase in NetCDF format.
 See ~abinit/scripts/post_processing/merge_ddb_nc.py
 for a merging utility.

COPYRIGHT

 Copyright (C) 1999-2018 ABINIT group (XG,MT,GA)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

  filename=name of the DDB.nc file to be written.
  Crystal <type(crystal_t)>=info on the crystal
  natom = number of atoms in the unit cell.
  mpert = maximum number of perturbations.
  qpt(3)=curret q-point, in reduced coordinates.
  d2matr(2,3,mpert,3,mpert)=second-order derivative of the total energy
  blkflg(3,mpert,3,mpert)=mask telling whether each element is computed (1) or not (0).

OUTPUT

  Only writing

PARENTS

      ddb_interpolate,respfn

CHILDREN

SOURCE

 38 #if defined HAVE_CONFIG_H
 39 #include "config.h"
 40 #endif
 41 
 42 #include "abi_common.h"
 43 
 44 subroutine outddbnc (filename, mpert, d2matr, blkflg, qpt, Crystal)
 45 
 46  use defs_basis
 47  use defs_datatypes
 48  use defs_abitypes
 49  use m_errors
 50  use m_xmpi
 51  use m_profiling_abi
 52  use m_nctk
 53  use m_crystal
 54 #ifdef HAVE_NETCDF
 55  use netcdf
 56 #endif
 57 
 58  use m_crystal_io, only : crystal_ncwrite
 59 
 60 !This section has been created automatically by the script Abilint (TD).
 61 !Do not modify the following lines by hand.
 62 #undef ABI_FUNC
 63 #define ABI_FUNC 'outddbnc'
 64 !End of the abilint section
 65 
 66  implicit none
 67 
 68 !Arguments -------------------------------
 69 !scalars
 70  character(len=*),intent(in) :: filename
 71  integer,intent(in) :: mpert
 72  integer,intent(in) :: blkflg(3,mpert,3,mpert)
 73  real(dp),intent(in) :: d2matr(2,3,mpert,3,mpert)
 74  real(dp),intent(in) :: qpt(3)
 75  type(crystal_t), intent(in) :: Crystal
 76 
 77 !Local variables -------------------------
 78 !scalars
 79  integer :: natom
 80  integer :: ncid, ncerr
 81  integer :: cplex, cart_dir, one_dim
 82  integer :: ipert1, ipert2, idir1, idir2
 83  integer,allocatable :: dynmat_mask(:,:,:,:)
 84  integer,allocatable :: born_effective_charge_tensor_mask(:,:,:)
 85  real(dp),allocatable :: dynmat(:,:,:,:,:)
 86  real(dp),allocatable :: born_effective_charge_tensor(:,:,:)
 87 
 88 ! *********************************************************************
 89 
 90 #ifdef HAVE_NETCDF
 91 
 92  natom = Crystal%natom
 93 
 94  ABI_ALLOCATE(dynmat, (2,3,natom,3,natom))
 95  ABI_ALLOCATE(dynmat_mask, (3,natom,3,natom))
 96  ABI_ALLOCATE(born_effective_charge_tensor, (3,natom,3))
 97  ABI_ALLOCATE(born_effective_charge_tensor_mask, (3,natom,3))
 98 
 99  ! Initialize NetCDF file.
100  NCF_CHECK(nctk_open_create(ncid, filename, xmpi_comm_self))
101 
102  ! ------------------------------
103  ! Construct local DDB
104  ! ------------------------------
105 
106  ! Construct the dynamical matrix
107  do ipert1=1,natom
108    do idir1=1,3
109      do ipert2=1,natom
110        do idir2=1,3
111          dynmat_mask(idir1,ipert1,idir2,ipert2) = blkflg(idir1,ipert1,idir2,ipert2)
112          if (blkflg(idir1,ipert1,idir2,ipert2)==1) then
113            dynmat(1,idir1,ipert1,idir2,ipert2) = d2matr(1,idir1,ipert1,idir2,ipert2)
114            dynmat(2,idir1,ipert1,idir2,ipert2) = d2matr(2,idir1,ipert1,idir2,ipert2)
115          else
116            dynmat(1,idir1,ipert1,idir2,ipert2) = zero
117            dynmat(2,idir1,ipert1,idir2,ipert2) = zero
118          end if
119        end do
120      end do
121    end do
122  end do 
123 
124  ! Construct the Born effective charge tensor
125  ipert2 = natom + 2
126  do ipert1=1,natom
127    do idir1=1,3
128      do idir2=1,3
129        born_effective_charge_tensor_mask(idir1,ipert1,idir2) = blkflg(idir1,ipert1,idir2,ipert2)
130        if (blkflg(idir1,ipert1,idir2,ipert2)==1) then
131             ! This is a real quantity
132          born_effective_charge_tensor(idir1,ipert1,idir2) = d2matr(1,idir1,ipert1,idir2,ipert2)
133        else
134          born_effective_charge_tensor(idir1,ipert1,idir2) = zero
135        end if
136      end do
137    end do
138  end do 
139 
140  ! TODO: also store the dielectric matrix
141 
142  ! ------------------------------
143  ! Write crystal info
144  ! ------------------------------
145  ncerr = crystal_ncwrite(Crystal, ncid)
146  NCF_CHECK(ncerr)
147 
148 
149  ! ------------------------------
150  ! Write DDB
151  ! ------------------------------
152 
153  ! Write the dimensions specified by ETSF
154  one_dim = 1
155  cplex = 2
156  cart_dir = 3
157 
158  ncerr = nctk_def_dims(ncid, [&
159 & nctkdim_t('current_one_dim', one_dim), &
160 & nctkdim_t('number_of_atoms', natom), &
161 & nctkdim_t('number_of_cartesian_directions', cart_dir), &
162 & nctkdim_t('number_of_perturbations', mpert), &
163 & nctkdim_t('cplex',cplex)], defmode=.True.)
164  NCF_CHECK(ncerr)
165 
166  ! Create the arrays
167  ncerr = nctk_def_arrays(ncid, [&
168  nctkarr_t('atomic_masses_amu', "dp", 'number_of_atom_species'),&
169  nctkarr_t('q_point_reduced_coord', "dp", 'number_of_cartesian_directions'),&
170  nctkarr_t('second_derivative_of_energy', "dp", 'cplex, &
171 & number_of_cartesian_directions, number_of_atoms, &
172 & number_of_cartesian_directions, number_of_atoms'), &
173  nctkarr_t('second_derivative_of_energy_mask', "i", '&
174 & number_of_cartesian_directions, number_of_atoms, &
175 & number_of_cartesian_directions, number_of_atoms'), &
176  nctkarr_t('born_effective_charge_tensor', "dp", '&
177 & number_of_cartesian_directions, number_of_atoms, &
178 & number_of_cartesian_directions'), &
179  nctkarr_t('born_effective_charge_tensor_mask', "i", ' &
180 & number_of_cartesian_directions, number_of_atoms, &
181 & number_of_cartesian_directions')])
182  NCF_CHECK(ncerr)
183 
184 ! Write data
185  NCF_CHECK(nctk_set_datamode(ncid))
186  NCF_CHECK(nf90_put_var(ncid, vid('atomic_masses_amu'), Crystal%amu))
187  NCF_CHECK(nf90_put_var(ncid, vid('q_point_reduced_coord'), qpt))
188  NCF_CHECK(nf90_put_var(ncid, vid('second_derivative_of_energy'), dynmat))
189  NCF_CHECK(nf90_put_var(ncid, vid('second_derivative_of_energy_mask'), dynmat_mask))
190  NCF_CHECK(nf90_put_var(ncid, vid('born_effective_charge_tensor'), born_effective_charge_tensor))
191  NCF_CHECK(nf90_put_var(ncid, vid('born_effective_charge_tensor_mask'), born_effective_charge_tensor_mask))
192 
193  ! Close file
194  NCF_CHECK(nf90_close(ncid))
195 
196  ! Deallocate stuff
197 
198  ABI_FREE(dynmat)
199  ABI_FREE(dynmat_mask)
200  ABI_FREE(born_effective_charge_tensor)
201  ABI_FREE(born_effective_charge_tensor_mask)
202 
203 #else
204  MSG_ERROR("NETCDF support required to write DDB.nc file.")
205 #endif
206 
207  contains
208    integer function vid(vname) 
209 
210 
211 !This section has been created automatically by the script Abilint (TD).
212 !Do not modify the following lines by hand.
213 #undef ABI_FUNC
214 #define ABI_FUNC 'vid'
215 !End of the abilint section
216 
217    character(len=*),intent(in) :: vname
218    vid = nctk_idname(ncid, vname)
219  end function vid
220 
221 end subroutine outddbnc