TABLE OF CONTENTS


ABINIT/wvl_projectors_free [ Functions ]

[ Top ] [ Functions ]

NAME

 wvl_projectors_free

FUNCTION

 Freeing routine.

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (DC)
 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

OUTPUT

SIDE EFFECTS

  proj <type(wvl_projectors_type)>=projectors informations in a wavelet basis.

PARENTS

      gstate,wvl_wfsinp_reformat

CHILDREN

      free_dft_psp_projectors

SOURCE

 32 #if defined HAVE_CONFIG_H
 33 #include "config.h"
 34 #endif
 35 
 36 #include "abi_common.h"
 37 
 38 subroutine wvl_projectors_free(proj)
 39 
 40  use defs_basis
 41  use defs_wvltypes
 42  use m_profiling_abi
 43  use m_errors
 44 #if defined HAVE_BIGDFT
 45   use BigDFT_API, only : free_DFT_PSP_projectors,deallocate_gwf
 46 #endif
 47 
 48 !This section has been created automatically by the script Abilint (TD).
 49 !Do not modify the following lines by hand.
 50 #undef ABI_FUNC
 51 #define ABI_FUNC 'wvl_projectors_free'
 52 !End of the abilint section
 53 
 54  implicit none
 55 
 56 !Arguments ------------------------------------
 57 !scalars
 58  type(wvl_projectors_type),intent(inout) :: proj
 59 
 60 !Local variables -------------------------
 61 #if defined HAVE_BIGDFT
 62  integer :: ii
 63 #endif
 64 
 65   ! *********************************************************************
 66 
 67 #if defined HAVE_BIGDFT
 68 
 69  call free_DFT_PSP_projectors(proj%nlpsp)
 70 
 71  if (allocated(proj%G)) then
 72    do ii=1,size(proj%G)
 73 !    MT dec 2014: cannot call bigdft deallocation routine
 74 !    because content of proj%G datastructure was created
 75 !    without f_malloc (without memory profiling).
 76 !    call deallocate_gwf(proj%G(ii))
 77      if (associated(proj%G(ii)%ndoc)) then
 78        ABI_DEALLOCATE(proj%G(ii)%ndoc)
 79      end if
 80      if (associated(proj%G(ii)%nam)) then
 81        ABI_DEALLOCATE(proj%G(ii)%nam)
 82      end if
 83      if (associated(proj%G(ii)%nshell)) then
 84        ABI_DEALLOCATE(proj%G(ii)%nshell)
 85      end if
 86      if (associated(proj%G(ii)%psiat)) then
 87        ABI_DEALLOCATE(proj%G(ii)%psiat)
 88      end if
 89      if (associated(proj%G(ii)%xp)) then
 90        ABI_DEALLOCATE(proj%G(ii)%xp)
 91      end if
 92    end do
 93    ABI_DATATYPE_DEALLOCATE(proj%G)
 94  end if
 95 
 96 #else
 97  BIGDFT_NOTENABLED_ERROR()
 98  if (.false.) write(std_out,*) proj%nlpsp
 99 #endif
100 
101 end subroutine wvl_projectors_free