TABLE OF CONTENTS
ABINIT/mblktyp5 [ Functions ]
NAME
mblktyp5
FUNCTION
This routine merges the derivative databases of type 5: second-order eigenvalue derivatives why is this separate from mblktyp1? Should be merged at some point for consistency
COPYRIGHT
Copyright (C) 1999-2018 ABINIT group (XG,MT,SP) This file is distributed under the terms of the GNU General Public Licence, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
NOTES
The heading of the database is read, then the heading of the temporary database to be added is read, the code check their compatibility, and create a new database that mixes the old and the temporary ones.
INPUTS
chkopt=option for consistency checks between DDB files codename=MRGDDB ddbun=define input and output unit numbers dscrpt=description of the output file filnam=name of input or output file mddb=maximum number of databases (cannot be made dynamic) nddb=number of input DDBs
OUTPUT
msym=maximum number of symmetry elements in space group Merge the file
PARENTS
mrgddb
CHILDREN
ddb_free,ddb_hdr_compare,ddb_hdr_free,ddb_hdr_open_read ddb_hdr_open_write,ddb_malloc,ddb_write_blok,read_blok8,wrtout
SOURCE
47 #if defined HAVE_CONFIG_H 48 #include "config.h" 49 #endif 50 51 #include "abi_common.h" 52 53 54 subroutine mblktyp5 (chkopt,ddbun,dscrpt,filnam,mddb,msym,nddb,vrsddb) 55 56 use defs_basis 57 use m_xmpi 58 use m_errors 59 use m_profiling_abi 60 use m_ddb 61 use m_ddb_hdr 62 63 !This section has been created automatically by the script Abilint (TD). 64 !Do not modify the following lines by hand. 65 #undef ABI_FUNC 66 #define ABI_FUNC 'mblktyp5' 67 use interfaces_14_hidewrite 68 !End of the abilint section 69 70 implicit none 71 72 !Arguments ------------------------------- 73 !scalars 74 integer,intent(in) :: ddbun,mddb,nddb,vrsddb 75 integer,intent(out) :: msym 76 character(len=fnlen),intent(in) :: dscrpt,filnam(mddb+1) 77 78 !Local variables ------------------------- 79 !scalars 80 !Define input and output unit numbers: 81 integer,parameter :: ddbuntmp=3 82 integer :: chkopt,choice,dimekb,iblok,iblok1,iblok2 83 integer :: iddb,ii,lmnmax,matom 84 integer :: mband,mblktyp,mblok,mkpt,mpert,msize,mtypat 85 integer :: nblok,nblokt 86 integer :: temp,tmerge,usepaw 87 integer,allocatable :: mgblok(:) 88 real(dp),parameter :: qtol=2.0d-8 89 real(dp) :: diff 90 type(ddb_type) :: ddb 91 type(ddb_hdr_type) :: ddb_hdr, ddb_hdr8 92 !arrays 93 real(dp),allocatable :: blkval2(:,:,:,:),kpnt(:,:,:) 94 character(len=500) :: message 95 96 ! ********************************************************************* 97 98 ! Make sure there is more than one ddb to be read 99 if(nddb==1)then 100 101 write(message, '(a,a,a,a,a)' )& 102 & 'The initialisation mode of MRGDDB, that uses nddb=1,',& 103 & 'has been disabled in version 2.1 of ABINIT.',& 104 & 'Action : you should use DDBs that include the symmetry',& 105 & 'information (and that can be used and merged without',& 106 & 'initialisation), or you should use ABINITv2.0.' 107 MSG_ERROR(message) 108 109 end if 110 111 !Evaluate the maximal dimensions of arrays 112 dimekb=0 ; matom=0 ; mband=0 ; mblok=0 ; mkpt=0 113 msize=0 ; mtypat=0 ; lmnmax=0 ; usepaw=0 ; mblktyp = 1 114 msym=192 115 116 do iddb=1,nddb 117 call ddb_hdr_open_read(ddb_hdr, filnam(iddb+1), ddbun, vrsddb,& 118 & dimonly=1) 119 120 mblok=mblok+ddb_hdr%nblok 121 mblktyp=max(mblktyp,ddb_hdr%mblktyp) 122 matom=max(matom,ddb_hdr%matom) 123 mkpt=max(mkpt,ddb_hdr%mkpt) 124 mtypat=max(mtypat,ddb_hdr%mtypat) 125 msym=max(msym,ddb_hdr%msym) 126 mband=max(mband,ddb_hdr%mband) 127 dimekb=max(dimekb,ddb_hdr%psps%dimekb) 128 lmnmax=max(lmnmax,ddb_hdr%psps%lmnmax) 129 usepaw=max(usepaw,ddb_hdr%usepaw) 130 131 132 call ddb_hdr_free(ddb_hdr) 133 134 end do 135 136 mpert=matom+6 137 msize=3*mpert*3*mpert 138 if(mblktyp==3)msize=msize*3*mpert 139 140 !write(std_out,*),'msize',msize,'mpert',mpert,'mblktyp',mblktyp 141 call ddb_malloc(ddb,msize,mblok,matom,mtypat) 142 143 !Allocate arrays 144 ABI_ALLOCATE(mgblok,(mblok)) 145 146 147 !********************************************************************** 148 149 !Read the first database 150 151 write(std_out,*)' read the input derivative database information' 152 call ddb_hdr_open_read(ddb_hdr, filnam(2), ddbun, vrsddb, & 153 & matom=matom,mtypat=mtypat,mband=mband,mkpt=mkpt,& 154 & msym=msym,dimekb=dimekb,lmnmax=lmnmax,usepaw=usepaw) 155 156 ABI_ALLOCATE(blkval2,(2,msize,ddb_hdr%nband(1),mkpt)) 157 ABI_ALLOCATE(kpnt,(3,mkpt,mblok)) 158 159 nblok = ddb_hdr%nblok 160 161 if(nblok>=1)then 162 ! Read the blocks from the input database. 163 write(message, '(a,i5,a)' ) ' read ',nblok,' blocks from the input DDB ' 164 call wrtout(std_out,message,'COLL') 165 choice=1 166 do iblok=1,nblok 167 call read_blok8(ddb,iblok,ddb_hdr%nband(1),mpert,& 168 & msize,ddb_hdr%nkpt,ddbun,blkval2(1,1,1,1),kpnt(1,1,iblok)) 169 ! Setup merged indicator 170 mgblok(iblok)=1 171 end do 172 else 173 call wrtout(std_out,' No bloks in the first ddb ','COLL') 174 end if 175 !Close the first ddb 176 close(ddbun) 177 178 !********************************************* 179 180 !In case of merging of DDBs, iterate the reading 181 do iddb=2,nddb 182 183 ! Open the corresponding input DDB, 184 ! and read the database file informations 185 write(message, '(a,a,i6)' )ch10,& 186 & ' read the input derivative database number',iddb 187 call wrtout(std_out,message,'COLL') 188 189 call ddb_hdr_open_read(ddb_hdr8, filnam(iddb+1), ddbun, vrsddb, & 190 & matom=matom,mtypat=mtypat,mband=mband,mkpt=mkpt,& 191 & msym=msym,dimekb=dimekb,lmnmax=lmnmax,usepaw=usepaw) 192 193 if (chkopt==1)then 194 ! Compare the current DDB and input DDB information. 195 ! In case of an inconsistency, halt the execution. 196 write(message, '(a)' )' compare the current and input DDB information' 197 call wrtout(std_out,message,'COLL') 198 199 call ddb_hdr_compare(ddb_hdr, ddb_hdr8) 200 201 else if(chkopt==0)then 202 ! No comparison between the current DDB and input DDB information. 203 write(message, '(a)' )' no comparison between the current and input DDB information' 204 call wrtout(std_out,message,'COLL') 205 write(message, '(a,a,a)' )& 206 & 'No comparison/check is performed for the current and input DDB information ',& 207 & 'because argument --nostrict was passed to the command line. ',& 208 & 'Use at your own risk !' 209 MSG_COMMENT(message) 210 end if 211 212 call wrtout(std_out,' Will try to merge this input DDB with the current one.','COLL') 213 214 ! First estimate of the total number of bloks, and error 215 ! message if too large 216 write(message, '(a,i5)' ) ' Current number of bloks =',nblok 217 call wrtout(std_out,message,'COLL') 218 write(message, '(a,i5,a)' )' Will read ',ddb_hdr8%nblok,' blocks from the input DDB ' 219 call wrtout(std_out,message,'COLL') 220 nblokt=nblok+ddb_hdr8%nblok 221 if(nblokt>mblok)then 222 write(message, '(a,i0,a,a,a,i0,a)' )& 223 & 'The expected number of blocks',nblokt,' is larger than',ch10,& 224 & 'the maximum number of blocks',mblok,'.' 225 MSG_ERROR(message) 226 end if 227 228 ! Read the bloks from the temporary database, and close it. 229 ! Also setup the merging indicator 230 choice=1 231 do iblok=nblok+1,nblokt 232 call read_blok8(ddb,iblok,ddb_hdr8%nband(1),mpert,& 233 & msize,ddb_hdr8%nkpt,ddbun,blkval2(1,1,1,1),kpnt(1,1,iblok)) 234 mgblok(iblok)=1 235 end do 236 close(ddbun) 237 238 nblok=nblokt 239 write(message, '(a,i5)' ) ' Now, current number of bloks =',nblok 240 call wrtout(std_out,message,'COLL') 241 242 ! In certain cases, the different DDB will have different information 243 ! on the pseudos (depending on fullinit) 244 ! Here, we copy the information of the last DDB file, 245 ! only to make the tests pass... 246 ddb_hdr%psps%indlmn(:,:,:) = ddb_hdr8%psps%indlmn(:,:,:) 247 ddb_hdr%psps%pspso(:) = ddb_hdr8%psps%pspso(:) 248 ddb_hdr%psps%ekb(:,:) = ddb_hdr8%psps%ekb(:,:) 249 250 call ddb_hdr_free(ddb_hdr8) 251 252 end do 253 254 call wrtout(std_out,' All DDBs have been read ','COLL') 255 256 !********************************************************* 257 258 !Check the equality of blocks, and eventually merge them 259 260 if(nblok>=1)then 261 call wrtout(std_out,' check the equality of blocks, and eventually merge ','COLL') 262 do iblok2=2,nblok 263 do iblok1=1,iblok2-1 264 ! Check the block type identity 265 if(ddb%typ(iblok1)==ddb%typ(iblok2))then 266 ! Check the wavevector identities 267 diff=abs(ddb%qpt(1,iblok1)-ddb%qpt(1,iblok2)) 268 diff=diff+abs(ddb%qpt(2,iblok1)-ddb%qpt(2,iblok2)) 269 diff=diff+abs(ddb%qpt(3,iblok1)-ddb%qpt(3,iblok2)) 270 if(abs(diff)<qtol)mgblok(iblok2)=0 271 end if 272 end do 273 end do 274 275 ! Count the final number of bloks 276 tmerge=0 277 do ii=1,nblok 278 if(mgblok(ii)==1)tmerge=tmerge+1 279 end do 280 temp = nblok-tmerge 281 nblok=tmerge 282 283 ! Summarize the merging phase 284 write(message, '(i6,a,i6,a)' )& 285 & temp,' blocks are merged; the new DDB will have ',nblok,' blocks.' 286 call wrtout(std_out,message,'COLL') 287 end if 288 289 !********************************************************************** 290 291 write(message, '(a,a)' )' open the output database, write the',' preliminary information ' 292 call wrtout(std_out,message,'COLL') 293 294 ddb_hdr%dscrpt = trim(dscrpt) 295 ddb_hdr%nblok = nblok !nblokt 296 ddb_hdr%mblktyp = mblktyp 297 298 call ddb_hdr_open_write(ddb_hdr, filnam(1), ddbun, fullinit=1) 299 300 if(nddb>1)then 301 302 ! Write the whole database 303 call wrtout(std_out,' write the DDB ','COLL') 304 ii = 1 !unit indicator of what will be merged 305 ! Create a temporary file to decrease memory need. 306 do iddb=1,nddb 307 call ddb_hdr_open_read(ddb_hdr8, filnam(iddb+1), ddbuntmp, vrsddb) 308 309 do iblok=1,ddb_hdr8%nblok 310 if(mgblok(ii)==1) then 311 call read_blok8(ddb,ii,ddb_hdr8%nband(1),mpert,& 312 & msize,ddb_hdr8%nkpt,ddbuntmp,blkval2(:,:,:,:),kpnt(:,:,ii)) 313 choice=2 314 call ddb_write_blok(ddb,ii,choice,ddb_hdr%nband(1),mpert,& 315 & msize,ddb_hdr8%nkpt,ddbun,blkval2(:,:,:,:),kpnt(:,:,ii)) 316 else 317 write(message, '(a,i4,a,i4,a)' )& 318 & ' Bloc number',iblok,' of DDB ',iddb,& 319 & ' was merged, so do not write it' 320 call wrtout(std_out,message,'COLL') 321 end if 322 ii = ii+1 323 end do 324 close(ddbuntmp) 325 call ddb_hdr_free(ddb_hdr8) 326 end do !iddb=1,nddb 327 328 ! Also write summary of bloks at the end 329 write(ddbun, '(/,a)' )' List of bloks and their characteristics ' 330 choice=3 331 do iblok=1,nblokt 332 if(mgblok(iblok)==1)then 333 call ddb_write_blok(ddb,iblok,choice,ddb_hdr%nband(1),mpert,& 334 & msize,ddb_hdr%nkpt,ddbun) 335 end if 336 end do 337 338 end if 339 340 close (ddbun) 341 342 !********************************************************************* 343 344 call ddb_hdr_free(ddb_hdr) 345 346 !Deallocate arrays 347 ABI_DEALLOCATE(mgblok) 348 ABI_DEALLOCATE(blkval2) 349 ABI_DEALLOCATE(kpnt) 350 351 call ddb_free(ddb) 352 353 end subroutine mblktyp5