TABLE OF CONTENTS


ABINIT/m_args_gs [ Modules ]

[ Top ] [ Modules ]

NAME

  m_args_gs

FUNCTION

  This module provides the definition of the args_gs_type
  used to tranfer some arguments to GS calculations,
  especially those depending on the image of the cell.

COPYRIGHT

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

SOURCE

18 #if defined HAVE_CONFIG_H
19 #include "config.h"
20 #endif
21 
22 #include "abi_common.h"
23 
24 MODULE m_args_gs
25 
26  use defs_basis
27  use m_abicore
28  use m_errors
29 
30  implicit none
31 
32  private

m_args_gs/args_gs_free [ Functions ]

[ Top ] [ m_args_gs ] [ Functions ]

NAME

  args_gs_free

FUNCTION

  Clean and destroy a args_gs datastructure

INPUTS

OUTPUT

SIDE EFFECTS

  args_gs(:)=<type(args_gs_type)>=args_gs datastructure

SOURCE

156 subroutine args_gs_free(args_gs)
157 
158 !Arguments ------------------------------------
159 !arrays
160  type(args_gs_type),intent(inout) :: args_gs
161 !Local variables-------------------------------
162 
163 !************************************************************************
164 
165  !@args_gs_type
166 
167  args_gs%amu         => null()
168  args_gs%mixalch     => null()
169  args_gs%dmatpawu    => null()
170  args_gs%upawu       => null()
171  args_gs%jpawu       => null()
172  args_gs%rprimd_orig => null()
173 
174 end subroutine args_gs_free

m_args_gs/args_gs_init [ Functions ]

[ Top ] [ m_args_gs ] [ Functions ]

NAME

  args_gs_init

FUNCTION

  Init a args_gs datastructure

INPUTS

  amu(:)= mass of each atom type
  cellcharge= cell charge
  mixalch(:)= mixing coefficients to generate alchemical pseudo atoms
  dmatpawu(:,:,:,:)= fixed occupation matrix for correlated orbitals (DFT+U or DMFT only)
  upawu(:) =value of U for the DFT+U or DMFT approach
  jpawu(:) =value of J for the DFT+U or DMFT approach

OUTPUT

SIDE EFFECTS

  args_gs=<type(args_gs_type)>=args_gs datastructure

SOURCE

114 subroutine args_gs_init(args_gs,amu,cellcharge,mixalch,dmatpawu,upawu,jpawu,rprimd_orig)
115 
116 !Arguments ------------------------------------
117  real(dp), intent(in) :: cellcharge
118 !arrays
119  real(dp),intent(in),target :: amu(:),dmatpawu(:,:,:,:),jpawu(:),mixalch(:,:),upawu(:)
120  real(dp),intent(in),target :: rprimd_orig(:,:)
121  type(args_gs_type),intent(inout) :: args_gs
122 !Local variables-------------------------------
123 
124 !************************************************************************
125 
126  !@args_gs_type
127  args_gs%cellcharge  = cellcharge
128  args_gs%amu         => amu
129  args_gs%mixalch     => mixalch
130  args_gs%dmatpawu    => dmatpawu
131  args_gs%upawu       => upawu
132  args_gs%jpawu       => jpawu
133  args_gs%rprimd_orig => rprimd_orig
134 
135 end subroutine args_gs_init

m_args_gs/args_gs_type [ Types ]

[ Top ] [ m_args_gs ] [ Types ]

NAME

 args_gs_type

FUNCTION

 This structured datatype contains some arguments of a GS calculation,
 especially the ones depending on the "images".

SOURCE

45  type, public :: args_gs_type
46 
47 ! WARNING : if you modify this datatype, please check whether there might be creation/destruction/copy routines,
48 ! declared in another part of ABINIT, that might need to take into account your modification.
49 
50 ! Real number
51   real(dp) :: cellcharge
52    ! cell charge input variable, for one dataset and one image
53 
54 ! Real (real(dp)) arrays
55 
56   real(dp), pointer :: amu(:)
57    ! amu(ntypat)
58    ! mass of each atom type
59 
60   real(dp), pointer :: mixalch(:,:)
61    ! mixalch(ntypat)
62    ! Mixing coefficients to generate alchemical pseudo atoms
63 
64   real(dp), pointer :: dmatpawu(:,:,:,:)
65    ! dmatpawu(2*lpawu+1,2*lpawu+1,nspinor*nsppol,natpawu)
66    ! Fixed occupation matrix for correlated orbitals (DFT+U or DMFT only)
67 
68   real(dp), pointer :: upawu(:)
69    ! upawu(ntypat)
70    ! Value of U for the DFT+U or DMFT approach
71 
72   real(dp), pointer :: jpawu(:)
73    ! jpawu(ntypat)
74    ! Value of J for the DFT+U or DMFT approach
75 
76   real(dp), pointer :: rprimd_orig(:,:)
77    ! rprimd_orig(3,3)
78    ! Original primitive vectors (usually the input variable)
79 
80  end type args_gs_type
81 
82 !public procedures.
83  public :: args_gs_init
84  public :: args_gs_free