TABLE OF CONTENTS


ABINIT/m_scup_dataset [ Modules ]

[ Top ] [ Modules ]

NAME

  m_scup_dataset

FUNCTION

  module with the type of the input variables for scale_up 
  when initialized this is a subtype of multibinit_dtset_type

COPYRIGHT

  Copyright (C) 2014-2024 ABINIT group (AM)
  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

17 #if defined HAVE_CONFIG_H
18 #include "config.h"
19 #endif
20 
21 #include "abi_common.h"
22 
23 module m_scup_dataset
24 
25  use defs_basis
26  use m_abicore
27  use m_errors
28 
29  use m_parser, only : intagm
30  use m_symtk,  only : matr3inv
31  use m_bz_mesh
32 
33  implicit none
34 
35  private
36 
37  public :: scup_dtset_type
38  public :: scup_dtset_init
39  public :: scup_dtset_free
40  public :: outvars_scup
41  public :: invars10scup
42  public :: scup_kpath_new 
43  public :: scup_kpath_print

m_scup_dataset/invars10scup [ Functions ]

[ Top ] [ m_scup_dataset ] [ Functions ]

NAME

 invars10scup

FUNCTION

 Open input file for the multibinit code, then reads or echoes the input information 
 for SCALE UP part that is linked with multibinit.

INPUTS

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

OUTPUT

 scup_dtset <type(scup_dtset_type)> = datatype with all the input variables

NOTES

 Should be executed by one processor only.

SOURCE

285 subroutine invars10scup(scup_dtset,lenstr,string)
286 
287 implicit none 
288 !Arguments -------------------------------
289 !scalars
290  integer,intent(in) :: lenstr
291  character(len=*),intent(in) :: string
292  type(scup_dtset_type),intent(inout) :: scup_dtset
293 
294 !Local variables -------------------------
295 !Dummy arguments for subroutine 'intagm' to parse input file
296 !Set routine version number here:
297 !scalars
298  integer :: ii,jdtset,marr,tread
299  character(len=500) :: message
300 !arrays
301  integer,allocatable :: intarr(:)
302  real(dp),allocatable :: dprarr(:)
303  !tmp integer to transfer to logicals 
304  integer :: tmp_int
305 
306 !*********************************************************************
307 
308  marr=30
309  ABI_MALLOC(intarr,(marr))
310  ABI_MALLOC(dprarr,(marr))
311 
312  jdtset=1
313 
314 !=====================================================================
315 ! Initialize Dataset with default values
316 !=====================================================================
317 
318 call scup_dtset_init(scup_dtset) 
319 tmp_int=0 
320 
321 !=====================================================================
322 !start reading in dimensions and non-dependent variables
323 !=====================================================================
324 
325 !A 
326 !B 
327 !C 
328 !D 
329 !E 
330  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_elec_model',tread,'INT')
331  if(tread==1) tmp_int=intarr(1)
332  if(tmp_int<0 .or. tmp_int>1 )then
333    write(message, '(a,I3,a,a,a,a,a)' )&
334 &   'scup_elec_model is',tmp_int,', but the only allowed values',ch10,&
335 &   'are 0 and 1.',ch10,&
336 &   'Action: correct scup_elec_model in your input file.'
337    ABI_ERROR(message)
338  end if
339  if(tmp_int == 1) scup_dtset%scup_elec_model = .TRUE.
340  tmp_int = 0 
341 
342 !F 
343  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_freezden',tread,'INT')
344  if(tread==1) tmp_int=intarr(1)
345  if(tmp_int<0 .or. tmp_int>1 )then
346    write(message, '(a,I3,a,a,a,a,a)' )&
347 &   'scup_freezden is',tmp_int,', but the only allowed values',ch10,&
348 &   'are 0 and 1.',ch10,&
349 &   'Action: correct scup_freezden in your input file.'
350    ABI_ERROR(message)
351  end if
352  if(tmp_int == 1) scup_dtset%scup_freezden = .TRUE.
353  tmp_int = 0 
354 
355  !G 
356 !H 
357 !I 
358  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_initorbocc',tread,'INT')
359  if(tread==1) tmp_int=intarr(1)
360  if(tmp_int<0 .or. tmp_int>1 )then
361    write(message, '(a,I3,a,a,a,a,a)' )&
362 &   'scup_initorbocc is',tmp_int,', but the only allowed values',ch10,&
363 &   'are 0 and 1.',ch10,&
364 &   'Action: correct scup_initorbocc in your input file.'
365    ABI_ERROR(message)
366  end if
367  if(tmp_int == 1) scup_dtset%scup_initorbocc = .TRUE.
368  tmp_int = 0 
369 
370  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_ismagnetic',tread,'INT')
371  if(tread==1) tmp_int=intarr(1)
372  if(tmp_int<0 .or. tmp_int>1 )then
373    write(message, '(a,I3,a,a,a,a,a)' )&
374 &   'scup_ismagnetic is',tmp_int,', but the only allowed values',ch10,&
375 &   'are 0 and 1.',ch10,&
376 &   'Action: correct scup_ismagnetic in your input file.'
377    ABI_ERROR(message)
378  end if
379  if(tmp_int == 1) scup_dtset%scup_ismagnetic = .TRUE.
380  tmp_int = 0 
381 
382  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_istddft',tread,'INT')
383  if(tread==1) tmp_int=intarr(1)
384  if(tmp_int<0 .or. tmp_int>1 )then
385    write(message, '(a,I3,a,a,a,a,a)' )&
386 &   'scup_istddft is',tmp_int,', but the only allowed values',ch10,&
387 &   'are 0 and 1.',ch10,&
388 &   'Action: correct scup_istddft in your input file.'
389    ABI_ERROR(message)
390  end if
391  if(tmp_int == 1) scup_dtset%scup_istddft = .TRUE.
392  tmp_int = 0 
393 
394 
395 !J
396 !K 
397  call intagm(dprarr,intarr,jdtset,marr,3,string(1:lenstr),'scup_ksamp',tread,'INT')
398  if(tread==1) scup_dtset%scup_ksamp(1:3)=intarr(1:3)
399  do ii=1,3
400    if(scup_dtset%scup_ksamp(ii)<1)then
401      write(message, '(a,i0,a,i0,a,a,a,i0,a)' )&
402 &     'scup_ksamp(',ii,') is',scup_dtset%scup_ksamp(ii),', which is lower than 1 .',ch10,&
403 &     'Action: correct scup_ksamp(',ii,') in your input file.'
404      ABI_ERROR(message)
405    end if
406  end do
407 
408 !L 
409 !M
410  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_maxscfstep',tread,'INT')
411  if(tread==1) scup_dtset%scup_maxscfstep=intarr(1)
412  if(scup_dtset%scup_maxscfstep<=0)then
413    write(message, '(a,I3,a,a,a,a,a)' )&
414 &   'scup_maxscfstep is',scup_dtset%scup_maxscfstep,', but the only allowed values',ch10,&
415 &   'greater than 0',ch10,&
416 &   'Action: correct scup_maxscfstep in your input file.'
417    ABI_ERROR(message)
418  end if
419 
420 !N
421 
422  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_ndivsm',tread,'INT')
423  if(tread==1) scup_dtset%scup_ndivsm=intarr(1)
424  if(scup_dtset%scup_ndivsm<0 )then
425    write(message, '(a,I3,a,a,a,a,a)' )&
426 &   'scup_ndivsm is',scup_dtset%scup_ndivsm,', but the only allowed values',ch10,&
427 &   'are positiv.',ch10,&
428 &   'Action: correct scup_ndivsm in your input file.'
429    ABI_ERROR(message)
430  end if
431 
432 
433  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_nspeck',tread,'INT')
434  if(tread==1) scup_dtset%scup_nspeck=intarr(1)
435  if(scup_dtset%scup_nspeck<0 )then
436    write(message, '(a,I3,a,a,a,a,a)' )&
437 &   'scup_nspeck is',scup_dtset%scup_nspeck,', but the only allowed values',ch10,&
438 &   'are positiv.',ch10,&
439 &   'Action: correct scup_nspeck in your input file.'
440    ABI_ERROR(message)
441  end if
442 
443 !O
444 !P 
445  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_printbands',tread,'INT')
446  if(tread==1) tmp_int=intarr(1)
447  if(tmp_int<0 .or. tmp_int>1 )then
448    write(message, '(a,I3,a,a,a,a,a)' )&
449 &   'scup_printbands is',tmp_int,', but the only allowed values',ch10,&
450 &   'are 0 and 1.',ch10,&
451 &   'Action: correct scup_printbands in your input file.'
452    ABI_ERROR(message)
453  end if
454  if(tmp_int == 1) scup_dtset%scup_printbands = .TRUE.
455  tmp_int = 0 
456 
457  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_printeigv',tread,'INT')
458  if(tread==1) tmp_int=intarr(1)
459  if(tmp_int<0 .or. tmp_int>1 )then
460    write(message, '(a,I3,a,a,a,a,a)' )&
461 &   'scup_printeigv is',tmp_int,', but the only allowed values',ch10,&
462 &   'are 0 and 1.',ch10,&
463 &   'Action: correct scup_printeigv in your input file.'
464    ABI_ERROR(message)
465  end if
466  if(tmp_int == 1) scup_dtset%scup_printeigv = .TRUE.
467  tmp_int = 0 
468 
469  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_printeltic',tread,'INT')
470  if(tread==1) tmp_int=intarr(1)
471  if(tmp_int<0 .or. tmp_int>1 )then
472    write(message, '(a,I3,a,a,a,a,a)' )&
473 &   'scup_printeltic is',tmp_int,', but the only allowed values',ch10,&
474 &   'are 0 and 1.',ch10,&
475 &   'Action: correct scup_printeltic in your input file.'
476    ABI_ERROR(message)
477  end if
478  if(tmp_int == 1) scup_dtset%scup_printeltic = .TRUE.
479  tmp_int = 0 
480 
481  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_printgeom',tread,'INT')
482  if(tread==1) tmp_int=intarr(1)
483  if(tmp_int<0 .or. tmp_int>1 )then
484    write(message, '(a,I3,a,a,a,a,a)' )&
485 &   'scup_printgeom is',tmp_int,', but the only allowed values',ch10,&
486 &   'are 0 and 1.',ch10,&
487 &   'Action: correct scup_printgeom in your input file.'
488    ABI_ERROR(message)
489  end if
490  if(tmp_int == 1) scup_dtset%scup_printgeom = .TRUE.
491  tmp_int = 0 
492 
493  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_printniter',tread,'INT')
494  if(tread==1) scup_dtset%scup_printniter=intarr(1)
495  if(scup_dtset%scup_printniter<0 )then
496    write(message, '(a,I3,a,a,a,a,a)' )&
497 &   'scup_printniter is',scup_dtset%scup_printniter,', but the only allowed values',ch10,&
498 &   'are positiv',ch10,&
499 &   'Action: correct scup_printniter in your input file.'
500    ABI_ERROR(message)
501  end if
502 
503  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_printorbocc',tread,'INT')
504  if(tread==1) tmp_int=intarr(1)
505  if(tmp_int<0 .or. tmp_int>1 )then
506    write(message, '(a,I3,a,a,a,a,a)' )&
507 &   'scup_printorbocc is',tmp_int,', but the only allowed values',ch10,&
508 &   'are 0 and 1.',ch10,&
509 &   'Action: correct scup_printorbocc in your input file.'
510    ABI_ERROR(message)
511  end if
512  if(tmp_int == 1) scup_dtset%scup_printorbocc = .TRUE.
513  tmp_int = 0 
514 
515 !Q
516 !R
517 !S
518 
519  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_startpulay',tread,'INT')
520  if(tread==1) scup_dtset%scup_startpulay=intarr(1)
521  if(scup_dtset%scup_startpulay<3)then
522    write(message, '(a,I3,a,a,a,a,a)' )&
523 &   'scup_startpulay is',scup_dtset%scup_startpulay,', but the only allowed values',ch10,&
524 &   'are greater than 3',ch10,&
525 &   'Action: correct scup_startpulay in your input file.'
526    ABI_ERROR(message)
527  end if
528 
529  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_scfmixing',tread,'DPR')
530  if(tread==1) scup_dtset%scup_scfmixing=dprarr(1)
531  if(scup_dtset%scup_scfmixing<0)then
532    write(message, '(a,f10.2,a,a,a,a,a)' )&
533 &   'scup_scfmixing is',scup_dtset%scup_scfmixing,', but the only allowed value',ch10,&
534 &   'is superior to 0.',ch10,&
535 &   'Action: correct scup_scfmixing in your input file.'
536    ABI_ERROR(message)
537  end if
538 
539  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_scfthresh',tread,'DPR')
540  if(tread==1) scup_dtset%scup_scfthresh=dprarr(1)
541  if(scup_dtset%scup_scfthresh <= 0)then
542    write(message, '(a,f10.2,a,a,a,a,a)' )&
543 &   'scup_scfthresh is',scup_dtset%scup_scfthresh,', but the only allowed value',ch10,&
544 &   'is superior to 0.',ch10,&
545 &   'Action: correct scup_scfthresh in your input file.'
546    ABI_ERROR(message)
547  end if
548 
549  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_smearing',tread,'DPR')
550  if(tread==1) scup_dtset%scup_smearing=dprarr(1)
551  if(scup_dtset%scup_smearing < 0)then
552    write(message, '(a,f10.2,a,a,a,a,a)' )&
553 &   'scup_smearing is',scup_dtset%scup_smearing,', but the only allowed value',ch10,&
554 &   'is superior to or equal to 0.',ch10,&
555 &   'Action: correct scup_smearing in your input file.'
556    ABI_ERROR(message)
557  end if
558 
559 !T
560  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'scup_tcharge',tread,'DPR')
561  if(tread==1) scup_dtset%scup_tcharge=dprarr(1)
562  if(scup_dtset%scup_tcharge<0)then
563    write(message, '(a,f10.2,a,a,a,a,a)' )&
564 &   'scup_tcharge is',scup_dtset%scup_tcharge,', but the only allowed value',ch10,&
565 &   'is superior to 0.',ch10,&
566 &   'Action: correct scup_tcharge in your input file.'
567    ABI_ERROR(message)
568  end if
569 
570 
571 !U 
572 !V 
573 !W
574 !X
575 !Y
576 !Z 
577 
578 !=====================================================================
579 !start reading dimension dependent variables
580 !=====================================================================
581 
582 
583 !A
584 !B
585 !C
586 !D
587 !E
588 !F
589 !G
590 !H
591 !I
592 !J
593 !K
594 !L
595 !M
596 !N
597 !O
598 !P
599 !Q
600 !R
601 !S
602    !Allocate
603    if(scup_dtset%scup_printbands)then 
604      ABI_MALLOC(scup_dtset%scup_speck,(3,scup_dtset%scup_nspeck))
605     
606    call intagm(dprarr,intarr,jdtset,marr,3*scup_dtset%scup_nspeck,string(1:lenstr),'scup_speck',tread,'DPR')
607      if(tread==1)then
608        scup_dtset%scup_speck(:,:)=reshape( dprarr(1:3*scup_dtset%scup_nspeck), [3,scup_dtset%scup_nspeck])
609      else
610        write(message,'(5a)') &
611 &       'When scup_printbands is asked, scup_speck must be initialized ',ch10,&
612 &       'in the input file, which is not the case.',ch10,&
613 &       'Action: initialize scup_speck in your input file, or change printbands.'
614        ABI_ERROR(message)
615      end if
616    end if
617 
618 
619 !T
620 !U 
621 !V 
622 !W
623 !X
624 !Y
625 !Z 
626 
627 !=======================================================================
628 !Finished reading in variables - deallocate
629 !=======================================================================
630 
631  ABI_FREE(dprarr)
632  ABI_FREE(intarr)
633 
634 end subroutine invars10scup

m_scup_dataset/outvars_scup [ Functions ]

[ Top ] [ m_scup_dataset ] [ Functions ]

NAME

 outvars_scup

FUNCTION

 Takes as an input the input dtset for scup and echoes it to 
 the output

INPUTS

 multibinit_dtset <type(multibinit_dtset_type)> datatype with all the input variables
 nunit=unit number for input or output

OUTPUT

  (only writing)

NOTES

 Should be executed by one processor only.

SOURCE

207 subroutine outvars_scup(scup_dtset,nunit)
208 
209  implicit none
210 
211 !Arguments -------------------------------
212 !scalars
213  integer,intent(in) :: nunit
214  type(scup_dtset_type),intent(in) :: scup_dtset
215 
216 !Local variables -------------------------
217 !Set routine version number here:
218 !scalars
219 !integers for printing
220  integer :: int_inorb=0,int_mgn=0,int_tddft=0,int_pband=0
221  integer :: int_peigv=0,int_peltic=0,int_pgeom=0,int_porbocc=0
222  integer :: int_freezden=0
223 !Character for defining fromat string 
224 !*********************************************************************
225    
226    !Check if logicals are true, if yes set integers to one for printing
227    if(scup_dtset%scup_initorbocc)   int_inorb    =1
228    if(scup_dtset%scup_ismagnetic)   int_mgn      =1
229    if(scup_dtset%scup_istddft)      int_tddft    =1
230    if(scup_dtset%scup_printbands)   int_pband    =1
231    if(scup_dtset%scup_printeigv)    int_peigv    =1
232    if(scup_dtset%scup_printeltic)   int_peltic   =1
233    if(scup_dtset%scup_printgeom)    int_pgeom    =1
234    if(scup_dtset%scup_printorbocc)  int_porbocc  =1
235    if(scup_dtset%scup_freezden)     int_freezden =1
236     
237    !Print  
238    write(nunit,'(a)')'Variables for SCALE-UP electronic model :'
239    write(nunit,'(1x,a16,3I3)')    '      scup_ksamp',scup_dtset%scup_ksamp
240    write(nunit,'(1x,a16,F7.3)')   '    scup_tcharge',scup_dtset%scup_tcharge
241    write(nunit,'(1x,a16,I3)')     ' scup_initorbocc',int_inorb  
242    write(nunit,'(1x,a16,I3)')     ' scup_ismagnetic',int_mgn    
243    write(nunit,'(1x,a16,I3)')     '    scup_istddft',int_tddft  
244    write(nunit,'(1x,a16,I3)')     ' scup_printbands',int_pband     
245    write(nunit,'(1x,a16,I3)')     '  scup_printeigv',int_peigv   
246    write(nunit,'(1x,a16,I3)')     ' scup_printeltic',int_peltic   
247    write(nunit,'(1x,a16,I3)')     '  scup_printgeom',int_pgeom    
248    write(nunit,'(1x,a16,I3)')     ' scup_printniter',scup_dtset%scup_printniter
249    write(nunit,'(1x,a16,I3)')     'scup_printorbocc',int_porbocc 
250    write(nunit,'(1x,a16,I3)')     '   scup_freezden',int_freezden
251    write(nunit,'(1x,a16,I3)')     '     scup_nspeck',scup_dtset%scup_nspeck
252    write(nunit,'(1x,a16,I3)')     '     scup_ndivsm',scup_dtset%scup_ndivsm
253    write(nunit,'(1x,a16,F7.3)')   '  scup_scfmixing',scup_dtset%scup_scfmixing 
254    write(nunit,'(1x,a16,ES10.2)') '  scup_scfthresh',scup_dtset%scup_scfthresh 
255    write(nunit,'(1x,a16,ES10.2)') '   scup_smearing',scup_dtset%scup_smearing 
256    write(nunit,'(1x,a16,I3)')     ' scup_startpulay',scup_dtset%scup_startpulay
257    write(nunit,'(1x,a16,I3)')     ' scup_maxscfstep',scup_dtset%scup_maxscfstep
258 
259 
260 end subroutine outvars_scup

m_scup_dataset/scup_dtset_free [ Functions ]

[ Top ] [ m_scup_dataset ] [ Functions ]

NAME

  scup_dtset_free

FUNCTION

  deallocate remaining arrays in the scup_dtset datastructure

INPUTS

  scup_dtset <type(scup_dtset_type)> = scup_dataset structure

 OUTPUTS
  scup_dtset <type(scup_dtset_type)> = scup_dataset structure

SOURCE

165 subroutine scup_dtset_free(scup_dtset)
166 
167  implicit none
168 
169 !Arguments ------------------------------------
170 !scalars
171  type(scup_dtset_type), intent(inout) :: scup_dtset
172 
173 ! *************************************************************************
174 
175 
176 if(allocated(scup_dtset%scup_speck))then 
177         ABI_FREE(scup_dtset%scup_speck)
178 endif
179 
180 call scup_dtset%scup_kpath%free
181 
182 
183 end subroutine scup_dtset_free

m_scup_dataset/scup_dtset_init [ Functions ]

[ Top ] [ m_scup_dataset ] [ Functions ]

NAME

 scup_dtset_init

FUNCTION

 Init the scup_dtset type

INPUTS

OUTPUT

 scup_dtset <type(scup_dtset_type)> = datatype with all the input variables

NOTES

 Should be executed by one processor only.

SOURCE

113 subroutine scup_dtset_init(scup_dtset)
114 
115 implicit none
116 
117 !Arguments -------------------------------
118 !scalars
119  type(scup_dtset_type),intent(inout) :: scup_dtset
120 !Local variables -------------------------
121 !scalars
122 !arrays
123 !-----------------------------------------
124 
125 scup_dtset%scup_ndivsm      =  0
126 scup_dtset%scup_nspeck      =  0 
127 scup_dtset%scup_elec_model  = .FALSE.
128 scup_dtset%scup_ksamp       =  (/ 1, 1, 1 /)
129 scup_dtset%scup_tcharge     =  0
130 scup_dtset%scup_initorbocc  = .FALSE. 
131 scup_dtset%scup_ismagnetic  = .FALSE. 
132 scup_dtset%scup_istddft     = .FALSE.
133 scup_dtset%scup_printbands  = .FALSE.  
134 scup_dtset%scup_printeigv   = .FALSE. 
135 scup_dtset%scup_printeltic  = .FALSE.  
136 scup_dtset%scup_printgeom   = .FALSE.  
137 scup_dtset%scup_printniter  =  0  
138 scup_dtset%scup_printorbocc = .FALSE. 
139 scup_dtset%scup_freezden    = .FALSE. 
140 scup_dtset%scup_scfmixing   =  0.3 
141 scup_dtset%scup_scfthresh   =  tol6
142 scup_dtset%scup_smearing    =  0.00091873313 ! Room Temperature in Hartree
143 scup_dtset%scup_startpulay  =  3
144 scup_dtset%scup_maxscfstep  =  100
145 
146 end subroutine  scup_dtset_init

m_scup_dataset/scup_dtset_type [ Types ]

[ Top ] [ m_scup_dataset ] [ Types ]

NAME

 scup_dtset_type

FUNCTION

 The scup_dtset_type structured datatype
 gathers all input variables and options for the SCALE UP part 
 that is linked with multibinit

SOURCE

57  type scup_dtset_type 
58 
59 !Integer 
60  integer :: scup_nspeck
61  integer :: scup_ndivsm
62  integer :: scup_printniter
63  integer :: scup_startpulay 
64  integer :: scup_maxscfstep
65 !Logicals 
66  logical :: scup_elec_model
67  logical :: scup_initorbocc
68  logical :: scup_ismagnetic 
69  logical :: scup_istddft    
70  logical :: scup_printbands
71  logical :: scup_printeigv
72  logical :: scup_printeltic 
73  logical :: scup_printgeom 
74  logical :: scup_printorbocc
75  logical(1) :: scup_freezden
76 !Real 
77  real*8   :: scup_tcharge 
78  real*8   :: scup_scfmixing
79  real*8   :: scup_scfthresh 
80  real*8   :: scup_smearing
81 !Integer Array
82  integer :: scup_ksamp(3)
83 
84 !Real Array 
85  real(dp),allocatable :: scup_speck(:,:)
86 
87 !Kpath Type 
88  type(kpath_t) :: scup_kpath
89 
90  end type scup_dtset_type

m_scup_dataset/scup_kpath_new [ Functions ]

[ Top ] [ m_scup_dataset ] [ Functions ]

NAME

 scup_kpath_init

FUNCTION

 Initialize the kpath and all other variables SCALE UP needs to plot electronic bands 
 along kpath

INPUTS

 speck = array with special k-points along the path 
 gprimd = reciprocal latice vectors of cell 
 ndivsm = number of divisions for smallest segment 

OUTPUT

 scup_kpath <type(kpath_t)> = kpath_t with all information about kpath

NOTES

 Should be executed by one processor only.

  m_bz_mesh/kpath_new

SOURCE

663 subroutine scup_kpath_new(speck,rprimd,ndivsm,scup_kpath)
664 
665  implicit none
666 !Arguments -------------------------------
667 !scalars
668  integer,intent(in) :: ndivsm
669 
670 !arrays 
671  real(dp),intent(in) :: speck(:,:),rprimd(3,3)
672  type(kpath_t), intent(out) :: scup_kpath
673 
674 !Local variables -------------------------
675 !Dummy arguments for subroutine 'intagm' to parse input file
676 !Set routine version number here:
677 !scalars
678  integer :: nspeck 
679 !arrays
680  integer,allocatable :: ndivs_tmp(:)
681  real(dp) :: gprimd(3,3)
682 
683 !*********************************************************************
684 
685 !Get gprimd 
686 call matr3inv(rprimd,gprimd)
687 
688 !Create Kpath 
689 scup_kpath = kpath_new(speck,gprimd,ndivsm)
690 
691 !Change size of scup_kpath%ndivs(:) variable 
692 !from nspeck-1 to nspeck and put 1 to first entry
693 nspeck = size(speck,2)
694 ABI_MALLOC(ndivs_tmp,(nspeck))
695 
696 !First entry is always 1
697  ndivs_tmp(1) = 1 
698 !Copy point/per segments
699  ndivs_tmp(2:) = scup_kpath%ndivs(:) 
700 
701 !Delete original segments 
702 ABI_FREE(scup_kpath%ndivs)
703 !Put new path 
704 ABI_MALLOC(scup_kpath%ndivs,(nspeck))
705 scup_kpath%ndivs = ndivs_tmp
706 !Free temporary array
707 ABI_FREE(ndivs_tmp)
708 
709 end subroutine scup_kpath_new

m_scup_dataset/scup_kpath_print [ Functions ]

[ Top ] [ m_scup_dataset ] [ Functions ]

NAME

 scup_kpath_print

FUNCTION

 Print info of kpath provide to SCALE UP 

INPUTS

 scup_kpath<tpye(kpath_t) = kpath_t with all information about kpath  

OUTPUT

 Only Printing

NOTES

 Should be executed by one processor only.

SOURCE

732 subroutine scup_kpath_print(scup_kpath)
733 
734  implicit none 
735 !Arguments ------------------------------------
736 !scalars
737 !arrays 
738  type(kpath_t), intent(in) :: scup_kpath 
739 !Local variables-------------------------------
740  integer :: unt
741 
742 ! *************************************************************************
743 
744  unt = std_out 
745  
746  write(unt,'(a)') ch10
747  write(unt,'(4a)') ' scup_printbands = 1. Printing of electronic bands active',ch10,&
748 &                  ' Kpath information below:',ch10
749 
750 call scup_kpath%print
751 
752 end subroutine scup_kpath_print