TABLE OF CONTENTS


ABINIT/mod_prc_memory [ Modules ]

[ Top ] [ Modules ]

NAME

 mod_prc_memory

FUNCTION

 This modules defines arrays and data used for the real-space kerker
 preconditionning of potential residuals.

COPYRIGHT

 Copyright (C) 2009-2018 ABINIT group (PMA).
 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 .

NOTES

  FIXME: this is highly non-kosher. Should be a datastructure which is declared dynamically
  MG: I completely agree! We don't use modules to share data and I don't see why we should
  break the rule here.

SOURCE

24 #if defined HAVE_CONFIG_H
25 #include "config.h"
26 #endif
27 
28 #include "abi_common.h"
29 
30 module mod_prc_memory
31 
32  use defs_basis
33  use m_abicore
34 
35  implicit none
36 
37 private
38 
39   real(dp),public,save,allocatable :: rdiemac(:)
40 
41   integer,save,public :: cycle=0 ! This is great! A global variable with the same name as a Fortran Statement!
42   real(dp),save,public :: energy_min
43 
44 public :: prc_mem_init
45 public :: prc_mem_free
46 
47 ! *********************************************************************
48 
49  contains

ABINIT/prc_mem_free [ Functions ]

[ Top ] [ Functions ]

NAME

 prc_mem_free

FUNCTION

 This subroutine deallocates the module's main component

PARENTS

      scfcv

CHILDREN

SOURCE

108 subroutine prc_mem_free()
109 
110 
111 !This section has been created automatically by the script Abilint (TD).
112 !Do not modify the following lines by hand.
113 #undef ABI_FUNC
114 #define ABI_FUNC 'prc_mem_free'
115 !End of the abilint section
116 
117 implicit none
118 
119 ! *********************************************************************
120 
121    if (allocated(rdiemac))  then
122      ABI_DEALLOCATE(rdiemac)
123    end if
124 
125  end subroutine prc_mem_free

ABINIT/prc_mem_init [ Functions ]

[ Top ] [ Functions ]

NAME

 prc_mem_init

FUNCTION

 This subroutine allocates the module's main component

PARENTS

      prcrskerker1

CHILDREN

SOURCE

66 subroutine prc_mem_init(nfft)
67 
68 
69 !This section has been created automatically by the script Abilint (TD).
70 !Do not modify the following lines by hand.
71 #undef ABI_FUNC
72 #define ABI_FUNC 'prc_mem_init'
73 !End of the abilint section
74 
75 implicit none
76 
77 !Arguments -------------------------------
78 integer, intent(in) :: nfft
79 !Local variables -------------------------
80 ! *********************************************************************
81 
82    if (.not. allocated(rdiemac))  then
83      ABI_ALLOCATE(rdiemac,(nfft))
84    end if
85    if(nfft.ne.size(rdiemac)) then ! This steps should be done over "istep" instead
86      ABI_DEALLOCATE(rdiemac)
87      ABI_ALLOCATE(rdiemac,(nfft))
88      cycle=0
89    end if
90 
91  end subroutine prc_mem_init