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