TABLE OF CONTENTS


ABINIT/m_band2eps_dataset [ Modules ]

[ Top ] [ Modules ]

NAME

  m_band2eps_dataset

FUNCTION

COPYRIGHT

  Copyright (C) 2014-2024 ABINIT group (XG,JCC,CL,MVeithen,XW,MJV)
  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

15 #if defined HAVE_CONFIG_H
16 #include "config.h"
17 #endif
18 
19 #include "abi_common.h"
20 
21 module m_band2eps_dataset
22 
23  use defs_basis
24  use m_abicore
25  use m_errors
26 
27  use m_parser,    only : intagm
28 
29  implicit none
30 
31  private
32 
33  public :: band2eps_dataset_type
34  public :: band2eps_dtset_free
35  public :: outvars_band2eps
36  public :: invars11

m_band2eps_dataset/band2eps_dataset_type [ Types ]

[ Top ] [ m_band2eps_dataset ] [ Types ]

NAME

 band2eps_dataset_type

FUNCTION

 The band2eps_dataset_type structured datatype
 gather all the input variables for the band2eps code.

SOURCE

51  type band2eps_dataset_type
52 
53 ! Variables should be declared on separated lines in order to reduce the occurence of git conflicts.
54 ! Since all these input variables are described in the band2eps_help.html
55 ! file, they are not described in length here ...
56 ! Integer
57   integer :: natom
58   integer :: cunit
59   integer :: nlines
60   integer :: ngrad
61   integer :: prtout
62   real(dp) :: min
63   real(dp) :: max
64 
65 
66   integer,allocatable :: nqline(:)
67 ! nqline(nlines)
68 
69   real(dp),allocatable :: scale(:)
70 ! scale(nlines)
71 
72   integer,allocatable :: red(:)
73 ! red(natom)
74 
75   integer,allocatable :: blue(:)
76 ! blue(natom)
77 
78   integer,allocatable :: green(:)
79 ! green(natom)
80 
81   character(len=6),allocatable :: qpoint_name(:)
82 ! qpoint_name(nlines+1)
83 
84 
85  end type band2eps_dataset_type

m_band2eps_dataset/band2eps_dtset_free [ Functions ]

[ Top ] [ m_band2eps_dataset ] [ Functions ]

NAME

   band2eps_dtset_free

FUNCTION

   deallocate remaining arrays in the band2eps_dtset datastructure

INPUTS

  band2eps_dtset = band2eps datastructure

NOTES

SOURCE

106 subroutine band2eps_dtset_free(band2eps_dtset)
107 
108 !Arguments ------------------------------------
109 !scalars
110  type(band2eps_dataset_type), intent(inout) :: band2eps_dtset
111 
112 ! *************************************************************************
113 
114  ABI_SFREE(band2eps_dtset%nqline)
115  ABI_SFREE(band2eps_dtset%scale)
116  ABI_SFREE(band2eps_dtset%red)
117  ABI_SFREE(band2eps_dtset%blue)
118  ABI_SFREE(band2eps_dtset%green)
119  ABI_SFREE(band2eps_dtset%qpoint_name)
120 
121 end subroutine band2eps_dtset_free

m_band2eps_dataset/invars11 [ Functions ]

[ Top ] [ m_band2eps_dataset ] [ Functions ]

NAME

 invars11

FUNCTION

 Open input file for the band2eps code, then reads or echoes the input information.

INPUTS

 lenstr=actual length of string
 natom=number of atoms, needed for atifc
 string*(*)=string of characters containing all input variables and data

OUTPUT

 band2eps_dtset= (derived datatype) contains all the input variables

NOTES

 Should be executed by one processor only.

SOURCE

148 subroutine invars11 (band2eps_dtset,lenstr,string)
149 
150 !Arguments -------------------------------
151 !scalars
152  integer,intent(in) :: lenstr
153  character(len=strlen),intent(in) :: string
154  type(band2eps_dataset_type),intent(inout) :: band2eps_dtset
155 
156 !Local variables -------------------------
157 !Dummy arguments for subroutine 'intagm' to parse input file
158 !Set routine version number here:
159 !scalars
160  integer :: ii,position,jdtset,marr,tread
161  character(len=500) :: message
162  character(len=fnlen) :: name_qpoint
163 !arrays
164  integer,allocatable :: intarr(:)
165  real(dp),allocatable :: dprarr(:)
166 
167 !*********************************************************************
168  marr=3
169  ABI_MALLOC(intarr,(marr))
170  ABI_MALLOC(dprarr,(marr))
171 
172  jdtset=1
173 
174 ! Need to get the number of atom
175  band2eps_dtset%natom = 0
176  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'natom',tread,'INT')
177  if(tread==1) band2eps_dtset%natom=intarr(1)
178  if (band2eps_dtset%natom <= tol10) then
179    write(message,'(a,I5,3a)' )&
180 &   'natom ',band2eps_dtset%natom,', is not allowed ',ch10,&
181 &   'Action: correct natom in your input file.'
182    ABI_ERROR(message)
183  end if
184 
185 ! Need to get the number of lines
186  band2eps_dtset%nlines = 0
187  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nlines',tread,'INT')
188  if(tread==1) band2eps_dtset%nlines=intarr(1)
189  if (band2eps_dtset%nlines <= tol10) then
190    write(message,'(a,I5,3a)' )&
191 &   'nlines ',band2eps_dtset%nlines,', is not allowed ',ch10,&
192 &   'Action: correct nlines in your input file.'
193    ABI_ERROR(message)
194  end if
195 
196 ! Need to get the number of lines
197  band2eps_dtset%cunit = 0
198  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'cunit',tread,'INT')
199  if(tread==1) band2eps_dtset%cunit=intarr(1)
200  if (band2eps_dtset%cunit < 1 .or. band2eps_dtset%cunit > 2) then
201    write(message,'(a,I5,3a)' )&
202 &   'cunit ',band2eps_dtset%cunit,', is not allowed ',ch10,&
203 &   'Action: correct cunit in your input file (1 or 2).'
204    ABI_ERROR(message)
205  end if
206 
207 ! Need to get the number of lines
208  band2eps_dtset%ngrad = 0
209  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ngrad',tread,'INT')
210  if(tread==1) band2eps_dtset%ngrad=intarr(1)
211  if (band2eps_dtset%ngrad < 0) then
212    write(message,'(a,I5,3a)' )&
213 &   'ngrad ',band2eps_dtset%ngrad,', is not allowed ',ch10,&
214 &   'Action: correct ngrad in your input file (positive value).'
215    ABI_ERROR(message)
216  end if
217 
218 ! Need to get the number of lines
219  band2eps_dtset%min = zero
220  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'min',tread,'DPR')
221  if(tread==1) band2eps_dtset%min=dprarr(1)
222 
223 ! Need to get the number of lines
224  band2eps_dtset%max = zero
225  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'max',tread,'DPR')
226  if(tread==1) band2eps_dtset%max=dprarr(1)
227 
228  ABI_MALLOC(band2eps_dtset%red,(band2eps_dtset%natom))
229  ABI_MALLOC(band2eps_dtset%blue,(band2eps_dtset%natom))
230  ABI_MALLOC(band2eps_dtset%green,(band2eps_dtset%natom))
231  band2eps_dtset%red(:) = 0
232  band2eps_dtset%blue(:) = 0
233  band2eps_dtset%green(:) = 0
234 
235 ! Read prtout
236  band2eps_dtset%prtout = 0
237  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'prtout',tread,'INT')
238  if(tread==1) band2eps_dtset%prtout=intarr(1)
239  if (band2eps_dtset%prtout < 0 .or. band2eps_dtset%prtout > 1) then
240    write(message,'(a,I5,3a)' )&
241 &   'prtout ',band2eps_dtset%prtout,', is not allowed ',ch10,&
242 &   'Action: correct prtout in your input file (0 or 1).'
243    ABI_ERROR(message)
244  end if
245 
246 !natom dimension
247  if(band2eps_dtset%natom > marr)then
248    marr = band2eps_dtset%natom
249    ABI_FREE(intarr)
250    ABI_FREE(dprarr)
251    ABI_MALLOC(intarr,(marr))
252    ABI_MALLOC(dprarr,(marr))
253  end if
254  call intagm(dprarr,intarr,jdtset,marr,band2eps_dtset%natom,string(1:lenstr),'red',tread,'INT')
255  if(tread==1) band2eps_dtset%red(1:band2eps_dtset%natom) = intarr(1:band2eps_dtset%natom)
256 
257  call intagm(dprarr,intarr,jdtset,marr,band2eps_dtset%natom,string(1:lenstr),'blue',tread,'INT')
258  if(tread==1) band2eps_dtset%blue(1:band2eps_dtset%natom) = intarr(1:band2eps_dtset%natom)
259 
260  call intagm(dprarr,intarr,jdtset,marr,band2eps_dtset%natom,string(1:lenstr),'green',tread,'INT')
261  if(tread==1) band2eps_dtset%green(1:band2eps_dtset%natom) = intarr(1:band2eps_dtset%natom)
262 
263 !nlines dimenstion
264  ABI_MALLOC(band2eps_dtset%nqline,(band2eps_dtset%nlines))
265  ABI_MALLOC(band2eps_dtset%scale,(band2eps_dtset%nlines))
266  band2eps_dtset%nqline(:) = 0
267  band2eps_dtset%scale(:) = zero
268 
269  if(band2eps_dtset%nlines > marr)then
270    marr = band2eps_dtset%nlines
271    ABI_FREE(intarr)
272    ABI_FREE(dprarr)
273    ABI_MALLOC(intarr,(marr))
274    ABI_MALLOC(dprarr,(marr))
275  end if
276 
277  call intagm(dprarr,intarr,jdtset,marr,band2eps_dtset%nlines,string(1:lenstr),'nqline',tread,'INT')
278  if(tread==1) band2eps_dtset%nqline(1:band2eps_dtset%nlines) = intarr(1:band2eps_dtset%nlines)
279 !REMOVE THE LAST LINE OF NQLINE
280   band2eps_dtset%nqline(band2eps_dtset%nlines) = band2eps_dtset%nqline(band2eps_dtset%nlines) - 1
281 
282  call intagm(dprarr,intarr,jdtset,marr,band2eps_dtset%nlines,string(1:lenstr),'scale',tread,'DPR')
283  if(tread==1) band2eps_dtset%scale(1:band2eps_dtset%nlines) = dprarr(1:band2eps_dtset%nlines)
284 
285 !nline+1 dimension
286  ABI_MALLOC(band2eps_dtset%qpoint_name,(band2eps_dtset%nlines+1))
287  band2eps_dtset%qpoint_name(:) = ""
288  position = index(string(1:lenstr),trim("QPOINT_NAME")) + 11
289  name_qpoint = trim(string(position:(position + &
290 &          len(band2eps_dtset%qpoint_name(1))*band2eps_dtset%nlines+1)))
291  read(name_qpoint,*) (band2eps_dtset%qpoint_name(ii),ii=1,band2eps_dtset%nlines+1)
292 
293 end subroutine invars11

m_band2eps_dataset/outvars_band2eps [ Functions ]

[ Top ] [ m_band2eps_dataset ] [ Functions ]

NAME

 outvars_band2eps

FUNCTION

 Open input file for the band2eps code, then
 echoes the input information.

INPUTS

 band2eps_dtset= (derived datatype) contains all the input variables
 nunit=unit number for input or output

OUTPUT

  (only writing)

NOTES

 Should be executed by one processor only.

SOURCE

319 subroutine outvars_band2eps (band2eps_dtset,nunit)
320 
321 !Arguments -------------------------------
322 !scalars
323  integer,intent(in) :: nunit
324  type(band2eps_dataset_type),intent(in) :: band2eps_dtset
325 
326 !Local variables -------------------------
327 !Set routine version number here:
328 !scalars
329  integer :: ii
330 
331 !*********************************************************************
332 
333 !Write the heading
334  write(nunit,'(80a,a)') ('=',ii=1,80),ch10
335  write(nunit, '(a,a)' )&
336 & ' -outvars_band2eps: echo values of input variables ----------------------',ch10
337 
338 !The flags
339  if(band2eps_dtset%natom/=0 )then
340    write(nunit,'(a)')' Informations :'
341    if(band2eps_dtset%natom/=0)write(nunit,'(3x,a9,3i10)')'  natom',band2eps_dtset%natom
342    if(band2eps_dtset%cunit/=0)write(nunit,'(3x,a9,3i10)')'  cunit',band2eps_dtset%cunit
343    if(band2eps_dtset%nlines/=0)write(nunit,'(3x,a9,3i10)')' nlines',band2eps_dtset%nlines
344    if(band2eps_dtset%min/=0)write(nunit,'(3x,a9,3f12.4)')'    Min',band2eps_dtset%min
345    if(band2eps_dtset%max/=0)write(nunit,'(3x,a9,3f12.4)')'    Max',band2eps_dtset%max
346    if(band2eps_dtset%ngrad/=0)write(nunit,'(3x,a9,3i10)')'  ngrad',band2eps_dtset%ngrad
347    write(nunit,'(3x,a9,6i10)')'    red',(band2eps_dtset%red(ii),ii=1,band2eps_dtset%natom)
348    write(nunit,'(3x,a9,6i10)')'   blue',(band2eps_dtset%blue(ii),ii=1,band2eps_dtset%natom)
349    write(nunit,'(3x,a9,6i10)')'  green',(band2eps_dtset%green(ii),ii=1,band2eps_dtset%natom)
350    write(nunit,'(3x,a9,8i10)')' nqline',(band2eps_dtset%nqline(ii),ii=1,band2eps_dtset%nlines)
351    write(nunit,'(3x,a9,8a)')'    point',(band2eps_dtset%qpoint_name(ii),ii=1,band2eps_dtset%nlines+1)
352  end if
353 
354  write(nunit,'(a,80a,a)') ch10,('=',ii=1,80),ch10
355 
356 end subroutine outvars_band2eps