00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112 #include "meam.h"
00113
00114
00115
00116 int MEAMFrame::readMEAM()
00117 {
00118 FILE *fp;
00119 char string[500], *p, *q;
00120
00121 LFile::SubHomeDir(potfile,potfile);
00122
00123 fp=fopen(potfile,"r");
00124 if(fp==NULL)
00125 {
00126 FATAL("MEAMFrame::readmeam file ("<<potfile<<") not found!");
00127 return -1;
00128 }
00129 fgets(string,500,fp);
00130 sscanf(string,"%lf",&zsmeam);
00131 INFO_Printf("zsmeam\t = %g\n",zsmeam);
00132
00133 fgets(string,500,fp);
00134 sscanf(string,"%lf",&alphas);
00135 INFO_Printf("alphas\t = %g\n",alphas);
00136
00137 fgets(string,500,fp);
00138 sscanf(string,"%lf %lf %lf %lf",betas,betas+1,betas+2,betas+3);
00139 INFO_Printf("betas\t = %g %g %g %g\n",betas[0],betas[1],betas[2],betas[3]);
00140
00141 fgets(string,500,fp);
00142 sscanf(string,"%lf %lf",&esubs,&asubs);
00143 INFO_Printf("esubs\t = %g asubs = %g\n",esubs,asubs);
00144
00145 fgets(string,500,fp);
00146 sscanf(string,"%lf %lf %lf %lf",ts,ts+1,ts+2,ts+3);
00147 INFO_Printf("ts\t = %g %g %g %g\n",ts[0],ts[1],ts[2],ts[3]);
00148
00149 fgets(string,500,fp);
00150 sscanf(string,"%lf",&rozros);
00151 INFO_Printf("rozros\t = %g\n",rozros);
00152
00153 fgets(string,500,fp);
00154 sscanf(string,"%d",&ibarr);
00155 INFO_Printf("ibarr\t = %d\n",ibarr);
00156
00157 fgets(string,500,fp);
00158 sscanf(string,"%lf",&rcutmeam);
00159 INFO_Printf("rcut\t = %g\n",rcutmeam);
00160
00161 fgets(string,500,fp);
00162 sscanf(string,"%lf %lf",&cmin,&cmax);
00163 INFO_Printf("cmin\t = %g cmax = %g\n",cmin,cmax);
00164
00165 fgets(string,500,fp);
00166 sscanf(string,"%lf %lf",&repuls,&attrac);
00167 INFO_Printf("repuls\t = %g attrac = %g\n",repuls,attrac);
00168
00169 fgets(string,500,fp);
00170 sscanf(string,"%lf",&legend);
00171 INFO_Printf("legend\t = %g (Legendre correction)\n",legend);
00172
00173 fgets(string,500,fp);
00174 sscanf(string,"%d",&noscr);
00175 INFO_Printf("noscr\t = %d\n",noscr);
00176
00177 fgets(string,500,fp);
00178 sscanf(string,"%lf",&alat);
00179 INFO_Printf("alat\t = %g\n",alat);
00180
00181 fgets(string,500,fp);
00182 sscanf(string,"%lf",&ielement);
00183 INFO_Printf("ielement\t = %g\n",ielement);
00184
00185 fgets(string,500,fp);
00186 p = strchr(string,'\'');
00187 if(p!=NULL)
00188 {
00189 q = strchr(p+1,'\'');
00190 if(q!=NULL)
00191 {
00192 q[0]=0;
00193 strcpy(lat,p+1);
00194 }
00195 }
00196 INFO_Printf("lat\t = %s\n",lat);
00197
00198 fgets(string,500,fp);
00199 p = strchr(string,'\'');
00200 if(p!=NULL)
00201 {
00202 q = strchr(p+1,'\'');
00203 if(q!=NULL)
00204 {
00205 q[0]=0;
00206 strcpy(elt,p+1);
00207 }
00208 }
00209 INFO_Printf("elt\t = %s\n",elt);
00210
00211 _RLIST = 1.1*rcutmeam;
00212 _SKIN = _RLIST - rcutmeam;
00213 rcutmeam2 = SQR(rcutmeam);
00214
00215 if(strcmp(lat,"dia")==0)
00216 {
00217 res = alat*sqrt(3.0)/4.0;
00218 INFO_Printf("diamond-cubic: res = %f\n",res);
00219 }
00220 else if (strcmp(lat,"fcc")==0)
00221 {
00222 res = alat/sqrt(2.0);
00223 INFO_Printf("FCC: res = %f\n",res);
00224 }
00225 else if (strcmp(lat,"bcc")==0)
00226 {
00227 res = alat*sqrt(3.0)/2.0;
00228 INFO_Printf("BCC: res = %f\n",res);
00229 }
00230 else
00231 {
00232 res = 0;
00233 FATAL("unknown lat ("<<lat<<"), res not set!");
00234 }
00235
00236 return 0;
00237 }
00238
00239 void MEAMFrame::rhoMEAM()
00240 {
00241
00242
00243
00244 int i, j, jpt, idx, jdx, itmp, jtmp, ktmp, ibar;
00245 Vector3 sij, rij, rcosv;
00246 double r2ij, rmagg, screenij, phifs, rej, rzj, drh;
00247 double betaj0, betaj1, betaj2, betaj3, rhj;
00248 double rogj1, rogj2, rogj3, rinve, rhocon;
00249 double t1, t2, t3, asub, esub, rozeroi, zmeam, z0meam;
00250 double rho0, rho1, rho2, rho3, rhobr, rhobar;
00251
00252
00253 rhocon = 1e-10;
00254
00255 for(i=0;i<_NP;i++)
00256 {
00257 atpe2b[i]=0;
00258 atpe3b[i]=0;
00259 _EPOT_IND[i]=0;
00260 rhotot[i]=0;
00261 embf [i]=0;
00262 c8a [i]=0;
00263 cg8c [i]=0;
00264
00265
00266 tav[i].clear();
00267 ag[i].clear();
00268 a8b[i].clear();
00269
00270
00271 b8c[i].clear();
00272
00273
00274 d8d[i].clear();
00275
00276 _F[i].clear();
00277 _VIRIAL_IND[i].clear();
00278 }
00279 _EPOT=0;
00280 _VIRIAL.clear();
00281
00282
00283
00284
00285
00286
00287 for(i=0;i<_NP;i++)
00288 {
00289 if(fixed[i]==-1) continue;
00290
00291 idx = species[i];
00292
00293 for(j=0;j<nn[i];j++)
00294 {
00295
00296 jpt=nindex[i][j];
00297 if(fixed[jpt]==-1) continue;
00298 jdx = species[jpt];
00299
00300 if(i==jpt) continue;
00301
00302 sij=_SR[jpt]-_SR[i];
00303 sij.subint();
00304 rij=_H*sij;
00305 r2ij=rij.norm2();
00306 rmagg=sqrt(r2ij);
00307
00308 #ifdef _EXTRACUTOFF
00309 if(rmagg>rcutmeam) continue;
00310 #endif
00311
00312 screenij = scrnab[i][j];
00313
00314 if(screenij<scnres) continue;
00315
00316
00317 phifs = phiid(rmagg,idx);
00318
00319
00320
00321
00322 atpe2b[i] += 0.5*phifs*screenij;
00323
00324
00325 rej = res;
00326 rzj = rozros;
00327 betaj0 = betas[0];
00328 betaj1 = betas[1];
00329 betaj2 = betas[2];
00330 betaj3 = betas[3];
00331
00332
00333 rhj = rhof(rmagg,betaj0,rej,rzj);
00334 rogj1 = rhof(rmagg,betaj1,rej,rzj);
00335 rogj2 = rhof(rmagg,betaj2,rej,rzj);
00336 rogj3 = rhof(rmagg,betaj3,rej,rzj);
00337
00338 rinve = 1./rmagg;
00339
00340
00341 rcosv = rij * rinve;
00342
00343 c8a[i] += rhj*screenij;
00344 for(itmp=0;itmp<3;itmp++)
00345 {
00346 t1 = screenij*rcosv[itmp];
00347 a8b[i][itmp] += rogj1*t1;
00348 ag [i][itmp] += rogj3*t1;
00349
00350 for(jtmp=0;jtmp<3;jtmp++)
00351 {
00352 t2 = t1*rcosv[jtmp];
00353 b8c[i][itmp][jtmp] += rogj2*t2;
00354
00355
00356 for(ktmp=0;ktmp<3;ktmp++)
00357 {
00358
00359 d8d[i].element[itmp][jtmp][ktmp] += rogj3*t2*rcosv[ktmp];
00360 }
00361 }
00362 }
00363
00364 tav[i][0] += ts[1]*rhj*screenij;
00365 tav[i][1] += ts[2]*rhj*screenij;
00366 tav[i][2] += ts[3]*rhj*screenij;
00367
00368 cg8c[i] += rogj2*screenij;
00369
00370 }
00371 }
00372
00373
00374 for(i=0;i<_NP;i++)
00375 {
00376 if(fixed[i]==-1) continue;
00377 rho0 = c8a[i];
00378
00379 if(rho0<=0)
00380 {
00381 embf[i] = 0;
00382 rhotot[i] = rhocon;
00383 rhsq[i].clear();
00384 }
00385 else
00386 {
00387 t1 = tav[i][0]/rho0;
00388 t2 = tav[i][1]/rho0;
00389 t3 = tav[i][2]/rho0;
00390 ibar = ibarr;
00391 asub = asubs;
00392 esub = esubs;
00393 rozeroi = rozros;
00394 zmeam = zsmeam;
00395
00396
00397 z0meam = zmeam;
00398
00399
00400
00401
00402 rho2 = -(cg8c[i]*cg8c[i])/3.0;
00403
00404 rho1 = rho3 = 0.0;
00405
00406
00407 for(itmp=0;itmp<3;itmp++)
00408 rho3 -= legend*SQR(ag[i][itmp]);
00409
00410 for(itmp=0;itmp<3;itmp++)
00411 {
00412 rho1 += SQR(a8b[i][itmp]);
00413 for(jtmp=0;jtmp<3;jtmp++)
00414 {
00415 rho2 += SQR(b8c[i][itmp][jtmp]);
00416 for(ktmp=0;ktmp<3;ktmp++)
00417 {
00418 rho3 += SQR(d8d[i].element[itmp][jtmp][ktmp]);
00419 }
00420 }
00421 }
00422
00423
00424 rhsq[i].set(rho1,rho2,rho3);
00425
00426
00427 drh = t1*rho1 + t2*rho2 + t3*rho3;
00428 rhobr = bar(rho0,drh,ibar,zmeam*rozeroi,dang1+i,dang2+i);
00429 rhobar= rhobr / (z0meam*rozeroi);
00430
00431 if(rhobr<0)
00432 {
00433 ERROR("rhobr("<<i<<") = "<<rhobr);
00434 }
00435 rhotot[i] = rhobr;
00436
00437 embf[i] = frhoi(rhobar,asub,esub);
00438
00439 atpe3b[i] = embf[i];
00440 }
00441 }
00442
00443 for(i=0;i<_NP;i++)
00444 {
00445 _EPOT_IND[i] = atpe2b[i] + atpe3b[i];
00446 _EPOT += _EPOT_IND[i];
00447 }
00448 }
00449
00450 void MEAMFrame::screen()
00451 {
00452 int i, j, k, jpt, kpt, idx, jdx, kdx, iid;
00453 Vector3 sij, rij, sik, rik, rjk;
00454 double r2ij, r2jk, r2ik, screenij;
00455 Vector3 *rij_table; double *r2ij_table;
00456
00457
00458
00459
00460
00461
00462
00463
00464
00465
00466 rij_table = (Vector3 *) malloc(sizeof(Vector3)*_NNM);
00467 r2ij_table = (double *) malloc(sizeof(double )*_NNM);
00468
00469 for(i=0;i<_NP;i++)
00470 {
00471 if(fixed[i]==-1) continue;
00472
00473 idx = species[i];
00474
00475
00476 for(j=0;j<nn[i];j++)
00477 {
00478 jpt=nindex[i][j];
00479 jdx = species[jpt];
00480 if(i==jpt) continue;
00481
00482 sij=_SR[jpt]-_SR[i];
00483 sij.subint();
00484 rij=_H*sij;
00485 r2ij=rij.norm2();
00486 rij_table[j] = rij;
00487 r2ij_table[j] = r2ij;
00488 }
00489
00490
00491
00492 for(j=0;j<nn[i];j++)
00493 {
00494
00495 jpt=nindex[i][j];
00496 if(fixed[jpt]==-1) continue;
00497 jdx = species[jpt];
00498
00499
00500 if(i>=jpt) continue;
00501
00502 iid=-1;
00503 for(k=0;k<nn[jpt];k++)
00504 {
00505 if(nindex[jpt][k]==i)
00506 {
00507 iid = k;
00508 break;
00509 }
00510 }
00511 if(iid==-1)
00512 FATAL("screen: i="<<i<<" j="<<j<<" jpt="<<jpt<<" iid="<<iid);
00513
00514
00515 scrnab[i][j] = 0;
00516 scrnab[jpt][iid] = 0;
00517
00518 rij = rij_table[j];
00519 r2ij = r2ij_table[j];
00520
00521 #ifdef _EXTRACUTOFF
00522
00523 if(r2ij>rcutmeam2) continue;
00524 #endif
00525
00526 screenij = rscrn(r2ij);
00527
00528
00529 if((noscr==1)||(screenij<=scnres)) continue;
00530 for(k=0;k<nn[i];k++)
00531 {
00532 if(k==j) continue;
00533 kpt = nindex[i][k];
00534 if(fixed[kpt]==-1) continue;
00535 kdx = species[kpt];
00536
00537 if(kpt==i) continue;
00538 if(kpt==jpt) continue;
00539
00540 rik = rij_table[k];
00541 r2ik = r2ij_table[k];
00542
00543 #ifdef _EXTRACUTOFF
00544
00545 if(r2ik>rcutmeam2) continue;
00546 #endif
00547
00548 rjk = rik - rij;
00549 r2jk=rjk.norm2();
00550
00551 #ifdef _EXTRACUTOFF
00552
00553 if(r2jk>rcutmeam2) continue;
00554 #endif
00555
00556
00557 screenij *= dscrn(r2ij,r2ik,r2jk,cmin,cmax);
00558
00559
00560
00561
00562
00563
00564 if(screenij<=scnres) break;
00565 }
00566 scrnab[i][j] = screenij;
00567 scrnab[jpt][iid] = screenij;
00568 }
00569 }
00570
00571 free(rij_table);
00572 free(r2ij_table);
00573 }
00574
00575 void MEAMFrame::dscreen()
00576 {
00577 int i, j, k, jpt, kpt, idx, jdx, kdx;
00578 Vector3 sij, rij, sik, rik, rjk, r2ij_p, r2ij_m, r2ik_p, r2ik_m;
00579 double r2ij, r2jk, r2ik, screenij;
00580 double hhmeam, h2meam, scrpx, scrpy, scrpz, scrmx, scrmy, scrmz;
00581
00582
00583
00584
00585
00586 if(noscr)
00587 {
00588 INFO("no need to call dscreen");
00589 return;
00590 }
00591
00592 hhmeam = SQR(hmeam);
00593 h2meam = 2.0 * hmeam;
00594
00595
00596 for(i=0;i<_NP;i++)
00597 {
00598 if(fixed[i]==-1) continue;
00599
00600 idx = species[i];
00601
00602 for(j=0;j<nn[i];j++)
00603 {
00604
00605 jpt=nindex[i][j];
00606 if(fixed[jpt]==-1) continue;
00607 jdx = species[jpt];
00608
00609 if(i==jpt) continue;
00610
00611 sij=_SR[jpt]-_SR[i];
00612 sij.subint();
00613 rij=_H*sij;
00614 r2ij=rij.norm2();
00615
00616 #ifdef _EXTRACUTOFF
00617
00618 if(r2ij>rcutmeam2) continue;
00619 #endif
00620
00621 screenij=scrnab[i][j];
00622
00623
00624
00625 if((screenij>1-scnres)||(screenij<scnres))
00626 {
00627 dscrnab[i][j*3] = dscrnab[i][j*3+1] = dscrnab[i][j*3+2] = 0;
00628 continue;
00629 }
00630
00631
00632 r2ij_p = rij*(-h2meam); r2ij_p += hhmeam + r2ij;
00633 r2ij_m = rij*( h2meam); r2ij_m += hhmeam + r2ij;
00634
00635 scrpx = rscrn(r2ij_p.x);
00636 scrpy = rscrn(r2ij_p.y);
00637 scrpz = rscrn(r2ij_p.z);
00638
00639 scrmx = rscrn(r2ij_m.x);
00640 scrmy = rscrn(r2ij_m.y);
00641 scrmz = rscrn(r2ij_m.z);
00642
00643
00644
00645
00646 for(k=0;k<nn[i];k++)
00647 {
00648 kpt = nindex[i][k];
00649 if(fixed[kpt]==-1) continue;
00650 kdx = species[kpt];
00651
00652 if(kpt==i) continue;
00653 if(kpt==jpt) continue;
00654
00655 sik=_SR[kpt]-_SR[i];
00656 sik.subint();
00657 rik=_H*sik;
00658 r2ik=rik.norm2();
00659
00660 #ifdef _EXTRACUTOFF
00661 if(r2ik>rcutmeam2) continue;
00662 #endif
00663 rjk = rik - rij;
00664 r2jk=rjk.norm2();
00665
00666 #ifdef _EXTRACUTOFF
00667 if(r2jk>rcutmeam2) continue;
00668 #endif
00669
00670
00671
00672 r2ik_p = rik*(-h2meam); r2ik_p += hhmeam + r2ik;
00673 r2ik_m = rik*( h2meam); r2ik_m += hhmeam + r2ik;
00674
00675 scrpx *= dscrn(r2ij_p.x,r2ik_p.x,r2jk,cmin,cmax);
00676 scrpy *= dscrn(r2ij_p.y,r2ik_p.y,r2jk,cmin,cmax);
00677 scrpz *= dscrn(r2ij_p.z,r2ik_p.z,r2jk,cmin,cmax);
00678
00679 scrmx *= dscrn(r2ij_m.x,r2ik_m.x,r2jk,cmin,cmax);
00680 scrmy *= dscrn(r2ij_m.y,r2ik_m.y,r2jk,cmin,cmax);
00681 scrmz *= dscrn(r2ij_m.z,r2ik_m.z,r2jk,cmin,cmax);
00682
00683
00684 }
00685 dscrnab[i][j*3+0] = (scrpx-scrmx)/(2*hmeam);
00686 dscrnab[i][j*3+1] = (scrpy-scrmy)/(2*hmeam);
00687 dscrnab[i][j*3+2] = (scrpz-scrmz)/(2*hmeam);
00688 }
00689 }
00690
00691 }
00692
00693 void MEAMFrame::dscrfor()
00694 {
00695 int i, j, k, jpt, kpt, idx, jdx, kdx, ibar, ib, itmp, jtmp, ktmp;
00696 Vector3 sij, rij, sik, rik, rjk, rkm, r2ij_p, r2ij_m, r2ik_p, r2ik_m, r2jk_p, r2jk_m;
00697 double r2ij, r2jk, r2ik, screenij;
00698 double hhmeam, h2meam;
00699 Vector3 rcosv, scrpk, scrmk, dscrnabk, df;
00700 double t1, t2, t3, s1, s2, s3, asub, esub, rozeroi, zmeam, z0meam, re, roz;
00701 double rho0, dfrho, drbd0, drbdAA, tmp0, rmagg, rinve;
00702 double rozero, beta0, beta1, beta2, beta3, rhs1, rhs2, rhs3, rhs4;
00703 double factor, ddd, phiij, scr1, drho0s, drho1s, drho1sg, drho2s, drho3s;
00704 double tmp1, tmp2, tmp3, tmp1g, term2, rml, rmln;
00705
00706
00707
00708
00709
00710
00711
00712
00713 if(noscr)
00714 {
00715 INFO("no need to call dscrfor");
00716 return;
00717 }
00718
00719 hhmeam = SQR(hmeam);
00720 h2meam = 2.0 * hmeam;
00721
00722 for(i=0;i<_NP;i++)
00723 {
00724 if(fixed[i]==-1) continue;
00725
00726 idx = species[i];
00727
00728 if(c8a[i]<=0)
00729 continue;
00730
00731
00732 t1 = tav[i][0]/c8a[i];
00733 t2 = tav[i][1]/c8a[i];
00734 t3 = tav[i][2]/c8a[i];
00735 ibar = ibarr;
00736 ib = abs(ibar);
00737 asub = asubs;
00738 esub = esubs;
00739 rozeroi = rozros;
00740 zmeam = zsmeam;
00741 z0meam = zmeam;
00742 s1 = 0.0;
00743 s2 = 0.0;
00744 s3 = 0.0;
00745
00746
00747 roz = rhotot[i] / (z0meam*rozeroi);
00748 rho0 = c8a[i];
00749 if(roz<0)
00750 {
00751 ERROR("roz["<<i<<"] = "<<roz);
00752 }
00753 dfrho = dfrhoi(roz,asub,esub)/(z0meam*rozeroi);
00754
00755
00756 drbd0 = rhotot[i]/rho0 + rho0*dang1[i];
00757
00758
00759 drbdAA = rho0*dang2[i];
00760
00761 tmp0 = (2.0/3.0) * cg8c[i];
00762
00763
00764
00765 for(j=0;j<nn[i];j++)
00766 {
00767
00768 jpt=nindex[i][j];
00769 if(fixed[jpt]==-1) continue;
00770 jdx = species[jpt];
00771
00772 if(i==jpt) continue;
00773
00774
00775 screenij=scrnab[i][j];
00776 if((screenij<=scnres)||(screenij>=1.0-scnres)) continue;
00777
00778 sij=_SR[jpt]-_SR[i];
00779 sij.subint();
00780 rij=_H*sij;
00781 r2ij=rij.norm2();
00782 rmagg=sqrt(r2ij);
00783
00784 #ifdef _EXTRACUTOFF
00785 if(rmagg>rcutmeam) continue;
00786 #endif
00787 rinve = 1./rmagg;
00788
00789
00790 rcosv = rij * rinve;
00791
00792 re = res;
00793 rozero = rozros;
00794 beta0 = betas[0];
00795 beta1 = betas[1];
00796 beta2 = betas[2];
00797 beta3 = betas[3];
00798
00799 rhs1 = rhof(rmagg,beta0,re,rozero);
00800 rhs2 = rhof(rmagg,beta1,re,rozero);
00801 rhs3 = rhof(rmagg,beta2,re,rozero);
00802 rhs4 = rhof(rmagg,beta3,re,rozero);
00803
00804 if(ibar<=0)
00805 {
00806 factor = 0.0;
00807 }
00808 else
00809 {
00810 factor = SQR(rho0/zmeam);
00811 }
00812
00813 s1 = s2 = s3 = 0.0;
00814
00815 ddd = (ts[1]-t1)*(rhsq[i][0]-s1*factor) +
00816 (ts[2]-t2)*(rhsq[i][1]-s2*factor) +
00817 (ts[3]-t3)*(rhsq[i][2]-s3*factor);
00818
00819 phiij = phiid(rmagg,idx);
00820
00821
00822
00823 for(k=0;k<nn[i];k++)
00824 {
00825 kpt = nindex[i][k];
00826 if(fixed[kpt]==-1) continue;
00827 kdx = species[kpt];
00828
00829 if(kpt==i) continue;
00830 if(kpt==jpt) continue;
00831
00832 sik=_SR[kpt]-_SR[i];
00833 sik.subint();
00834 rik=_H*sik;
00835 r2ik=rik.norm2();
00836
00837 #ifdef _EXTRACUTOFF
00838 if(r2ik>rcutmeam2) continue;
00839 #endif
00840 rjk = rik - rij;
00841 r2jk=rjk.norm2();
00842
00843 #ifdef _EXTRACUTOFF
00844 if(r2jk>rcutmeam2) continue;
00845 #endif
00846
00847 rkm = rik*2.0 - rij;
00848
00849
00850 scr1 = dscrn(r2ij,r2ik,r2jk,cmin,cmax);
00851
00852
00853 if((scr1<scnres)||(scr1>1.0-scnres)) continue;
00854
00855
00856 r2ik_p = rik*( h2meam); r2ik_p += r2ik + hhmeam;
00857 r2ik_m = rik*(-h2meam); r2ik_m += r2ik + hhmeam;
00858
00859 r2jk_p = rjk*( h2meam); r2jk_p += r2jk + hhmeam;
00860 r2jk_m = rjk*(-h2meam); r2jk_m += r2jk + hhmeam;
00861
00862 scrpk.x = dscrn(r2ij,r2ik_p.x,r2jk_p.x,cmin,cmax);
00863 scrpk.y = dscrn(r2ij,r2ik_p.y,r2jk_p.y,cmin,cmax);
00864 scrpk.z = dscrn(r2ij,r2ik_p.z,r2jk_p.z,cmin,cmax);
00865
00866 scrmk.x = dscrn(r2ij,r2ik_m.x,r2jk_m.x,cmin,cmax);
00867 scrmk.y = dscrn(r2ij,r2ik_m.y,r2jk_m.y,cmin,cmax);
00868 scrmk.z = dscrn(r2ij,r2ik_m.z,r2jk_m.z,cmin,cmax);
00869
00870
00871
00872 dscrnabk = scrpk - scrmk;
00873 dscrnabk *= (screenij/(scr1*hmeam*2));
00874
00875
00876 drho0s = rhs1;
00877 drho1s = 0.0;
00878 drho1sg = 0.0;
00879 drho2s = -rhs3 * tmp0;
00880 drho3s = 0.0;
00881
00882 for(itmp=0;itmp<3;itmp++)
00883 {
00884 tmp1 = 2.0 * a8b[i][itmp];
00885 tmp1g = 2.0 * ag[i][itmp];
00886 drho1s += tmp1*rhs2*rcosv[itmp];
00887 drho1sg += tmp1g*rhs4*rcosv[itmp];
00888
00889 for(jtmp=0;jtmp<3;jtmp++)
00890 {
00891 rml = rcosv[itmp]*rcosv[jtmp];
00892 tmp2 = 2.0 * b8c[i][itmp][jtmp];
00893 drho2s += tmp2 * rhs3 * rml;
00894
00895 for(ktmp=0;ktmp<3;ktmp++)
00896 {
00897 rmln = rml * rcosv[ktmp];
00898 tmp3 = 2.0 * d8d[i].element[itmp][jtmp][ktmp];
00899 drho3s += tmp3 * rhs4 * rmln;;
00900 }
00901 }
00902 }
00903
00904 term2 = t1*drho1s + t2*drho2s
00905 + t3*(drho3s-legend*drho1sg) + ddd*drho0s/rho0;
00906
00907
00908 df = dscrnabk * (0.5*phiij + dfrho * (drbd0*drho0s + drbdAA*term2));
00909
00910 _F[kpt] -= df;
00911
00912 _VIRIAL.addnvv(-0.5,df,rkm);
00913 _VIRIAL_IND[kpt].addnvv(-0.5,df,rkm);
00914 }
00915 }
00916 }
00917
00918 }
00919
00920
00921 void MEAMFrame::kraMEAM()
00922 {
00923 int i, j, ia, jpt, idx, jdx, ibar, ib, kk, itmp, jtmp, ktmp;
00924 double t1, t2, t3, asub, esub, rozeroi, zmeam, z0meam, s1, s2, s3, roz, rho0;
00925 double dfrho, drbd0, drbdAA, screenij, r2ij, rmagg, rinve;
00926 double beta0, beta1, beta2, beta3, re, rozero, rhs1, rhs2, rhs3, rhs4, drhs1, drhs2, drhs3, drhs4;
00927 double factor, ddd, tmp0;
00928 double drho0, drho1, drho1g, drho1sg, drho2, drho3;
00929 double drho0s, drho1s, drho2s, drho3s, tmp1, tmp2, tmp3, tmp4, tmp1g, term1, term2;
00930 double t1j, t2j, t3j, rozj, rho0j, dfrhoj, drbd0j, drbdAAj;
00931 double drho0j, drho1j, drho1gj, drho1sgj, drho2j, drho3j, dddj, tmp0j;
00932 double drho0sj, drho1sj, drho2sj, drho3sj, tmp1j, tmp2j, tmp3j, tmp4j, tmp1gj, term1j, term2j;
00933 double rml, rmln, drmllm, df, phip, phim, dphif, phiij;
00934 Vector3 sij, rij, rcosv, dfi, dfj, dfip, dfjp;
00935 Matrix33 drcos;
00936
00937
00938
00939
00940 for(i=0;i<_NP;i++)
00941 _F[i].clear();
00942
00943 for(i=0;i<_NP;i++)
00944 {
00945 if(fixed[i]==-1) continue;
00946
00947 idx = species[i];
00948
00949
00950 t1 = tav[i][0]/c8a[i];
00951 t2 = tav[i][1]/c8a[i];
00952 t3 = tav[i][2]/c8a[i];
00953 ibar = ibarr;
00954 ib = abs(ibar);
00955 asub = asubs;
00956 esub = esubs;
00957 rozeroi = rozros;
00958 zmeam = zsmeam;
00959 z0meam = zmeam;
00960 s1 = 0.0;
00961 s2 = 0.0;
00962 s3 = 0.0;
00963
00964
00965 roz = rhotot[i] / (z0meam*rozeroi);
00966 rho0 = c8a[i];
00967 if(roz<0)
00968 {
00969 ERROR("roz["<<i<<"] = "<<roz);
00970 }
00971 dfrho = dfrhoi(roz,asub,esub)/(z0meam*rozeroi);
00972
00973
00974 drbd0 = rhotot[i]/rho0 + rho0*dang1[i];
00975
00976
00977 drbdAA = rho0*dang2[i];
00978
00979
00980 for(j=0;j<nn[i];j++)
00981 {
00982
00983 jpt=nindex[i][j];
00984 if(fixed[jpt]==-1) continue;
00985 jdx = species[jpt];
00986
00987
00988 if(i>=jpt) continue;
00989
00990 screenij = scrnab[i][j];
00991 if(screenij<=scnres) continue;
00992
00993 sij=_SR[jpt]-_SR[i];
00994 sij.subint();
00995 rij=_H*sij;
00996
00997 r2ij=rij.norm2();
00998 rmagg=sqrt(r2ij);
00999
01000 #ifdef _EXTRACUTOFF
01001 if(rmagg>rcutmeam) continue;
01002 #endif
01003
01004 rinve = 1./rmagg;
01005
01006
01007 rcosv = rij * rinve;
01008
01009
01010 drcos = dcmij(rij,rmagg);
01011
01012
01013
01014 t1j = tav[jpt][0]/c8a[jpt];
01015 t2j = tav[jpt][1]/c8a[jpt];
01016 t3j = tav[jpt][2]/c8a[jpt];
01017
01018
01019 rozj = rhotot[jpt] / (z0meam*rozeroi);
01020 rho0j = c8a[jpt];
01021 if(rozj<0)
01022 {
01023 ERROR("roz["<<jpt<<"] = "<<rozj);
01024 }
01025 dfrhoj = dfrhoi(rozj,asub,esub)/(z0meam*rozeroi);
01026
01027
01028 drbd0j = rhotot[jpt]/rho0j + rho0j*dang1[jpt];
01029
01030
01031 drbdAAj = rho0j*dang2[jpt];
01032
01033 re = res;
01034 rozero = rozros;
01035 beta0 = betas[0];
01036 beta1 = betas[1];
01037 beta2 = betas[2];
01038 beta3 = betas[3];
01039
01040 rhs1 = rhof(rmagg,beta0,re,rozero);
01041 rhs2 = rhof(rmagg,beta1,re,rozero);
01042 rhs3 = rhof(rmagg,beta2,re,rozero);
01043 rhs4 = rhof(rmagg,beta3,re,rozero);
01044
01045 drhs1 = drhof(beta0,re,rhs1);
01046 drhs2 = drhof(beta1,re,rhs2);
01047 drhs3 = drhof(beta2,re,rhs3);
01048 drhs4 = drhof(beta3,re,rhs4);
01049
01050 if(ibar<=0)
01051 factor = 0.0;
01052 else
01053 factor = SQR(rho0/zmeam);
01054
01055 ddd = (ts[1]-t1)*(rhsq[i][0]-s1*factor) +
01056 (ts[2]-t2)*(rhsq[i][1]-s2*factor) +
01057 (ts[3]-t3)*(rhsq[i][2]-s3*factor);
01058
01059 tmp0 = (2.0/3.0) * cg8c[i];
01060
01061 dddj = (ts[1]-t1j)*(rhsq[jpt][0]-s1*factor) +
01062 (ts[2]-t2j)*(rhsq[jpt][1]-s2*factor) +
01063 (ts[3]-t3j)*(rhsq[jpt][2]-s3*factor);
01064
01065 tmp0j = (2.0/3.0) * cg8c[jpt];
01066
01067 for(kk=0;kk<3;kk++)
01068 {
01069 drho0 = drhs1*screenij*rcosv[kk];
01070 drho1 = 0.0;
01071 drho1g = 0.0;
01072 drho1sg = 0.0;
01073 drho2 = -tmp0*drhs3*screenij*rcosv[kk];
01074 drho3 = 0.0;
01075 drho0s = rhs1;
01076 drho1s = 0.0;
01077 drho2s = -rhs3*tmp0;
01078 drho3s = 0.0;
01079
01080 drho0j = drhs1*screenij*rcosv[kk];
01081 drho1j = 0.0;
01082 drho1gj = 0.0;
01083 drho1sgj = 0.0;
01084 drho2j = -tmp0j*drhs3*screenij*rcosv[kk];
01085 drho3j = 0.0;
01086 drho0sj = rhs1;
01087 drho1sj = 0.0;
01088 drho2sj = -rhs3*tmp0j;
01089 drho3sj = 0.0;
01090
01091 for(itmp=0;itmp<3;itmp++)
01092 {
01093 tmp1 = 2.0*a8b[i][itmp];
01094 tmp1g = 2.0*ag[i][itmp];
01095
01096 drho1 += tmp1*(drhs2*rcosv[kk]*rcosv[itmp] +
01097 rhs2*drcos[kk][itmp])*screenij;
01098 drho1g += tmp1g*(drhs4*rcosv[kk]*rcosv[itmp] +
01099 rhs4*drcos[kk][itmp])*screenij;
01100 drho1s += tmp1*rhs2*rcosv[itmp];
01101 drho1sg += tmp1g*rhs4*rcosv[itmp];
01102
01103 tmp1j = 2.0*a8b[jpt][itmp];
01104 tmp1gj = 2.0*ag[jpt][itmp];
01105
01106 drho1j += tmp1j*(drhs2*rcosv[kk]*rcosv[itmp] +
01107 rhs2*drcos[kk][itmp])*screenij;
01108 drho1gj += tmp1gj*(drhs4*rcosv[kk]*rcosv[itmp] +
01109 rhs4*drcos[kk][itmp])*screenij;
01110 drho1sj += tmp1j*rhs2*rcosv[itmp];
01111 drho1sgj += tmp1gj*rhs4*rcosv[itmp];
01112
01113 for(jtmp=0;jtmp<3;jtmp++)
01114 {
01115 rml = rcosv[itmp]*rcosv[jtmp];
01116 drmllm = drcos[kk][itmp]*rcosv[jtmp]+
01117 rcosv[itmp]*drcos[kk][jtmp];
01118 tmp2 = 2.0*b8c[i][itmp][jtmp];
01119 drho2 += tmp2*(drhs3*rcosv[kk]*rml + rhs3*drmllm)*screenij;
01120 drho2s += tmp2*rhs3*rml;
01121
01122 tmp2j = 2.0*b8c[jpt][itmp][jtmp];
01123 drho2j += tmp2j*(drhs3*rcosv[kk]*rml + rhs3*drmllm)*screenij;
01124 drho2sj += tmp2j*rhs3*rml;
01125
01126 for(ktmp=0;ktmp<3;ktmp++)
01127 {
01128 rmln = rml*rcosv[ktmp];
01129 tmp3 = 2.0*d8d[i].element[itmp][jtmp][ktmp];
01130 drho3 += tmp3*(drhs4*rcosv[kk]*rmln
01131 + rhs4*(drmllm*rcosv[ktmp] + rml*drcos[kk][ktmp]))
01132 * screenij;
01133 drho3s += tmp3*rhs4*rmln;
01134
01135 tmp3j = 2.0*d8d[jpt].element[itmp][jtmp][ktmp];
01136 drho3j += tmp3j*(drhs4*rcosv[kk]*rmln
01137 + rhs4*(drmllm*rcosv[ktmp] + rml*drcos[kk][ktmp]))
01138 * screenij;
01139 drho3sj += tmp3j*rhs4*rmln;
01140 }
01141 }
01142 }
01143
01144 term1 = t1*drho1 + t2*drho2 + t3*(drho3-legend*drho1g) + ddd*drho0/rho0;
01145
01146 term2 = t1*drho1s + t2*drho2s + t3*(drho3s-legend*drho1sg) + ddd*drho0s/rho0;
01147
01148 term1j = t1*drho1j - t2*drho2j + t3*(drho3j-legend*drho1gj) + dddj*drho0j/rho0j;
01149
01150 term2j = t1*drho1sj - t2*drho2sj + t3*(drho3sj-legend*drho1sgj) + dddj*drho0sj/rho0j;
01151
01152
01153 for(ia=0;ia<nn[jpt];ia++)
01154 if(nindex[jpt][ia]==i) break;
01155 if(ia==nn[jpt])
01156 FATAL(jpt<<"in neighbor list of "<<i<<" but the reverse it not true");
01157
01158 tmp4 = dscrnab[i][3*j+kk];
01159 df = dfrho*((drho0-drho0s*tmp4)*drbd0 + (term1-term2*tmp4)*drbdAA);
01160 dfi[kk]=df;
01161
01162 tmp4 = dscrnab[jpt][3*ia+kk];
01163 df = dfrho*((drho0+drho0s*tmp4)*drbd0 + (term1+term2*tmp4)*drbdAA);
01164 dfj[kk]=df;
01165
01166 tmp4j = dscrnab[jpt][3*ia+kk];
01167 df = dfrhoj*((-drho0j-drho0sj*tmp4j)*drbd0j + (term1j+term2j*tmp4j)*drbdAAj);
01168 dfip[kk]=df;
01169
01170 tmp4j = dscrnab[i][3*j+kk];
01171 df = dfrhoj*((-drho0j+drho0sj*tmp4j)*drbd0j + (term1j-term2j*tmp4j)*drbdAAj);
01172 dfjp[kk]=df;
01173 }
01174
01175
01176
01177
01178
01179
01180
01181
01182
01183
01184
01185 _F[i] += dfi; _F[jpt] -= dfj;
01186 _F[i] -= dfjp; _F[jpt] += dfip;
01187
01188 _VIRIAL.addnvv(-0.5,dfi,rij);
01189 _VIRIAL_IND[i].addnvv(-0.5,dfi,rij);
01190
01191 _VIRIAL.addnvv(-0.5,dfj,rij);
01192 _VIRIAL_IND[jpt].addnvv(-0.5,dfj,rij);
01193
01194 _VIRIAL.addnvv( 0.5,dfjp,rij);
01195 _VIRIAL_IND[i].addnvv( 0.5,dfjp,rij);
01196
01197 _VIRIAL.addnvv( 0.5,dfip,rij);
01198 _VIRIAL_IND[jpt].addnvv( 0.5,dfip,rij);
01199
01200
01201 phip = phiid(rmagg+hmeam,idx);
01202 phim = phiid(rmagg-hmeam,idx);
01203
01204 dphif=(phip-phim)/(2.*hmeam);
01205
01206 phiij=phiid(rmagg,idx);
01207
01208 dfi = rcosv*(dphif*screenij);
01209 dfi[0] -= phiij*dscrnab[i][3*j+0];
01210 dfi[1] -= phiij*dscrnab[i][3*j+1];
01211 dfi[2] -= phiij*dscrnab[i][3*j+2];
01212
01213
01214 _F[i] += dfi;
01215 _VIRIAL.addnvv(-0.5,dfi,rij);
01216 _VIRIAL_IND[i].addnvv(-0.5,dfi,rij);
01217
01218
01219
01220
01221 dfi = rcosv*(dphif*screenij*(-1.0));
01222 dfi[0] -= phiij*dscrnab[jpt][3*ia+0];
01223 dfi[1] -= phiij*dscrnab[jpt][3*ia+1];
01224 dfi[2] -= phiij*dscrnab[jpt][3*ia+2];
01225 _F[jpt] += dfi;
01226 _VIRIAL.addnvv( 0.5,dfi,rij);
01227 _VIRIAL_IND[jpt].addnvv( 0.5,dfi,rij);
01228 }
01229 }
01230 }
01231
01232
01233 double MEAMFrame::phiid(double r, int it)
01234 {
01235 double charge, re, alpha, xmjc, phizbl, frac, fdimer;
01236 double rozero, esub, repul, attra;
01237 double asub, ztmp, t1, t2, t3, z0meam, rh0, roz, zdimer;
01238 double beta0, beta1, beta2, beta3, rh1, rh2, rh3, arg, dum1, dum2;
01239 double phi;
01240 int ibar;
01241
01242
01243
01244 charge = ielement;
01245 re = res;
01246 alpha = alphas;
01247
01248 beta0=betas[0];
01249 beta1=betas[1];
01250 beta2=betas[2];
01251 beta3=betas[3];
01252 t1 =ts[1];
01253 t2 =ts[2];
01254 t3 =ts[3];
01255
01256 rozero=rozros;
01257 esub =esubs;
01258 repul =repuls;
01259 attra =attrac;
01260 asub =asubs;
01261 ibar =ibarr;
01262 ztmp =zsmeam;
01263
01264 z0meam=zbar(ibar,ztmp,lat,t1,t2,t3);
01265
01266
01267 xmjc=alpha*(r/re-1.0);
01268
01269
01270
01271
01272
01273
01274
01275
01276
01277
01278
01279
01280
01281
01282
01283
01284
01285
01286
01287
01288
01289
01290 rh0=rhof(r,beta0,re,1.0);
01291 rh1=rhof(r,beta1,re,1.0);
01292 rh2=rhof(r,beta2,re,1.0);
01293 rh3=rhof(r,beta3,re,1.0);
01294 if(xmjc<=xzbl0)
01295 {
01296 arg=t1*SQR(rh1)+(2.0/3.0)*t2*SQR(rh2)+t3*(1-legend)*SQR(rh3);
01297 zdimer=1;
01298 roz=bar(zdimer*rh0,arg,ibar,zdimer,&dum1,&dum2)/z0meam;
01299 fdimer=frhoi(roz,asub,esub);
01300
01301 phizbl=zbl(r,charge,charge);
01302
01303 if(xmjc<xzbl)
01304 frac=0.0;
01305 else
01306 frac=xcut((xmjc-xzbl0)/(xzbl-xzbl0));
01307 }
01308 else
01309 {
01310 frac=1;
01311 phizbl=0;
01312 fdimer=0;
01313 }
01314
01315
01316
01317
01318
01319
01320
01321
01322
01323
01324
01325
01326
01327
01328
01329
01330
01331
01332
01333
01334
01335 rh0=rhof(r,beta0,re,1.0);
01336
01337
01338 if(strcmp(lat,"dia")==0)
01339 {
01340 rh3=rhof(r,beta3,re,1.0);
01341
01342 arg=(32.0/9.0)*t3*SQR(rh3);
01343 rh0=bar(z0meam*rh0,arg,ibar,z0meam,&dum1,&dum2);
01344 roz=rh0/z0meam;
01345
01346 }
01347 else if(strcmp(lat,"hcp")==0)
01348 {
01349 rh3=rhof(r,beta3,re,1.0);
01350 arg=(1.0/3.0)*ts[3]*SQR(rh3);
01351 rh0=bar(z0meam*rh0,arg,ibar,z0meam,&dum1,&dum2);
01352 roz=rh0/z0meam;
01353 }
01354 else if(strcmp(lat,"dim")==0)
01355 {
01356 rh1=rhof(r,beta1,re,1.0);
01357 rh2=rhof(r,beta2,re,1.0);
01358 rh3=rhof(r,beta3,re,1.0);
01359
01360 arg=t1*SQR(rh1)+(2.0/3.0)*t2*SQR(rh2)+t3*(1-legend)*SQR(rh3);
01361 rh0=bar(z0meam*rh0,arg,ibar,z0meam,&dum1,&dum2);
01362 roz=rh0/z0meam;
01363 }
01364 else
01365 {
01366 FATAL("lat = "<<lat<<" not recognized");
01367 }
01368
01369
01370 phi=(2.0/ztmp)*( erose(r,re,alpha,esub,repul,attra)
01371 - frhoi(roz,asub,esub) );
01372
01373
01374
01375
01376
01377
01378
01379
01380
01381
01382 if(!enable_zbl_fdimer) fdimer=0.0;
01383
01384
01385
01386 phi=frac*phi+(1.0-frac)*(phizbl-2*fdimer);
01387
01388
01389
01390
01391
01392
01393
01394
01395
01396 return phi;
01397 }
01398
01399 double MEAMFrame::zbl(double r,double z1val,double z2val)
01400 {
01401 int i; double z1, z2;
01402 double c[4]={0.028171,0.28022,0.50986,0.18175};
01403 double d[4]={0.20162,0.40290,0.94229,3.1998};
01404 double azero = 0.4685;
01405 double cc = 14.3997;
01406 double x, a, ret;
01407
01408
01409
01410 z1 = round(z1val);
01411 z2 = round(z2val);
01412
01413 a=azero/(pow(z1,0.23)+pow(z2,0.23));
01414 ret=0.0;
01415
01416 x = r/a;
01417 for(i=0;i<4;i++)
01418 ret += c[i]*exp(-d[i]*x);
01419
01420 ret *= (z1*z2/r*cc);
01421
01422 return ret;
01423 }
01424
01425 double MEAMFrame::bar(double rho0,double A,int ibar,double z,double *dang1,double *dang2)
01426 {
01427 int ib;
01428 double ang, ex, x;
01429 double x0=-0.99; int n=99;
01430
01431 if(rho0<=0) return 1e-20;
01432
01433 ib = abs(ibar);
01434 if((ib==0)||(ib==4))
01435 {
01436 ang = 1+A/SQR(rho0);
01437
01438
01439
01440
01441
01442
01443
01444 if(ang<1+x0)
01445 {
01446 x=ang-1;
01447 ang=sqrt((1+x0)*pow(x0/x,n));
01448 *dang1=ang*n/rho0;
01449 *dang2=-ang*n/(2*A);
01450 }
01451 else
01452 {
01453 ang=sqrt(ang);
01454 *dang1 = -A/(rho0*rho0*rho0*ang);
01455 *dang2 = 1/(2*rho0*rho0*ang);
01456 }
01457 }
01458 else if(ib==1)
01459 {
01460 ang = exp(0.5*A/SQR(rho0));
01461 *dang1 = -ang*A/(rho0*rho0*rho0);
01462 *dang2 = ang/(2*rho0*rho0);
01463 }
01464 else if(ib==2)
01465 {
01466 ang = exp(0.5*A/SQR(z));
01467 *dang1 = 0;
01468 *dang2 = ang/(2.0*SQR(z));
01469 }
01470 else if(ib==3)
01471 {
01472 ex = exp(-A/SQR(rho0));
01473 ang = 2.0/(1.0 + ex);
01474 *dang1 = -SQR(ang)*ex*A/(rho0*rho0*rho0);
01475 *dang2 = SQR(ang)/2.0*ex/SQR(rho0);
01476 }
01477 return rho0*ang;
01478 }
01479
01480 inline double MEAMFrame::erose(double r,double re,double alpha,double esub,double repuls,double attrac)
01481 {
01482 double x, an3, E;
01483
01484
01485
01486
01487 x=alpha*(r/re-1.0);
01488 if(x>=0)
01489 an3=attrac;
01490 else
01491 an3=repuls;
01492
01493 E=-esub*(1.+x+an3*x*x*x/(r/re))*exp(-x);
01494
01495
01496
01497
01498
01499
01500
01501
01502
01503
01504
01505
01506
01507
01508
01509 return E;
01510 }
01511
01512
01513
01514
01515
01516
01517 Matrix33 MEAMFrame::dcmij(Vector3 rr,double rs)
01518 {
01519 int i, j; double r1, r13, tmp1;
01520 Matrix33 drcos;
01521
01522 r1 = 1.0 / rs;
01523 r13 = r1 * r1 * r1;
01524 for(i=0;i<3;i++)
01525 {
01526 tmp1 = r13 * rr[i];
01527 for(j=0;j<3;j++)
01528 {
01529 drcos[j][i] = -tmp1 * rr[j];
01530 }
01531 drcos[i][i] += r1;
01532 }
01533 return drcos;
01534 }
01535
01536 inline double MEAMFrame::zbar(int ibar,double z,char *lattice,double t1,double t2,double t3)
01537 {
01538 if(ibar<=0)
01539 return z;
01540 else
01541 {
01542 FATAL("ibar = "<<ibar<<" not implemented!");
01543 return 0;
01544 }
01545 }
01546
01547 inline double MEAMFrame::rhof(double r,double abc,double re,double rozero)
01548 {
01549
01550
01551
01552 return rozero*exp(-abc*(r/re-1.0));
01553 }
01554
01555 inline double MEAMFrame::frhoi(double rhotp,double asub,double esub)
01556 {
01557 if(rhotp>0)
01558 return asub*esub*rhotp*log(rhotp);
01559 else
01560 return asub*esub*rhotp*(-100.0);
01561 }
01562
01563 inline double MEAMFrame::dfrhoi(double rho,double asub,double esub)
01564 {
01565 if(rho>0)
01566 return asub*esub*(log(rho)+1.0);
01567 else
01568 return asub*esub*(-100.0);
01569 }
01570
01571 inline double MEAMFrame::rscrn(double r2ij)
01572 {
01573 double fac, rij, frcut=0.9;
01574
01575 if(enable_square_rscrn)
01576 {
01577
01578 fac=(r2ij-frcut*rcutmeam2)/((1.0-frcut)*rcutmeam2);
01579 }
01580 else
01581 {
01582
01583 rij = sqrt(r2ij);
01584 fac = 1.0-(rcutmeam-rij)/(1.0-frcut);
01585 }
01586 if(fac <= 0.0)
01587 return 1.0;
01588 else if(fac >= 1.0)
01589 return 0.0;
01590 else
01591 return xcut(fac);
01592 }
01593
01594 double MEAMFrame::dscrn(double rij, double rik, double rjk, double cmin, double cmax)
01595 {
01596 double dotil, dotjl, fac;
01597 double argb, argt, arg;
01598
01599 if((rik>sconst*rij)||(rjk>sconst*rij)) return 1.0;
01600
01601
01602 dotil = rik+rij-rjk;
01603 if(dotil<=2.0e-3) return 1.0;
01604 dotjl = rij+rjk-rik;
01605 if(dotjl<=2.0e-3) return 1.0;
01606 argb=SQR(rij)-SQR(rik-rjk);
01607 argt=4.0*(rik*rij)-SQR(rik-rjk+rij);
01608 arg=argt/argb;
01609
01610 if(arg>=cmax) return 1.0;
01611 else if(arg<=cmin) return 0.0;
01612 else
01613 {
01614 fac = (cmax-arg)/(cmax-cmin);
01615
01616 return xcut(fac);
01617 }
01618
01619 #if 0
01620 argb = 1.0 - SQR((rik-rjk)/rij);
01621 if(argb==0) return 1.0;
01622
01623 argt = 4.0*(rik/rij)-SQR((rik-rjk+rij)/rij);
01624 arg=argt/argb;
01625
01626
01627
01628
01629
01630
01631
01632
01633
01634
01635
01636
01637
01638
01639
01640
01641
01642
01643 if((arg>cmax)||(arg<0.0))
01644 return 1.0;
01645 else if(arg<cmin)
01646 return 0.0;
01647 else
01648 return xcut( (cmax-arg)/(cmax-cmin) );
01649 #endif
01650 }
01651
01652
01653
01654
01655
01656 void MEAMFrame::Alloc()
01657 {
01658 MDPARALLELFrame::Alloc();
01659 int i, size;
01660 size=_NP*allocmultiple;
01661
01662
01663 Realloc(atpe2b,double,size);
01664 Realloc(atpe3b,double,size);
01665 Realloc(rhotot,double,size);
01666 Realloc(embf,double,size);
01667 Realloc(c8a,double,size);
01668 Realloc(dang1,double,size);
01669 Realloc(dang2,double,size);
01670 Realloc(cg8c,double,size);
01671 Realloc(tav,Vector3,size);
01672 Realloc(rhsq,Vector3,size);
01673 Realloc(a8b,Vector3,size);
01674 Realloc(ag,Vector3,size);
01675 Realloc(b8c,Matrix33,size);
01676 Realloc(d8d,Matrix333,size);
01677
01678
01679 Realloc(scrnab_mem,char,size*_NNM*sizeof(double)+size*sizeof(double *));
01680 scrnab=(double **)(scrnab_mem+size*_NNM*sizeof(double));
01681 for(i=0;i<size;i++) scrnab[i]=(double *)(scrnab_mem+i*_NNM*sizeof(double));
01682
01683 Realloc(dscrnab_mem,char,size*_NNM*3*sizeof(double)+size*sizeof(double *));
01684 dscrnab=(double **)(dscrnab_mem+size*_NNM*3*sizeof(double));
01685 for(i=0;i<size;i++) dscrnab[i]=(double *)(dscrnab_mem+i*_NNM*3*sizeof(double));
01686
01687
01688 bindvar("atpe2b",atpe2b,DOUBLE);
01689 bindvar("atpe3b",atpe3b,DOUBLE);
01690 bindvar("rhotot",rhotot,DOUBLE);
01691 bindvar("embf",embf,DOUBLE);
01692 bindvar("c8a",c8a,DOUBLE);
01693 bindvar("dang1",dang1,DOUBLE);
01694 bindvar("dang2",dang2,DOUBLE);
01695 bindvar("cg8c",cg8c,DOUBLE);
01696 bindvar("tav",tav,DOUBLE);
01697 bindvar("rhsq",rhsq,DOUBLE);
01698 bindvar("a8b",a8b,DOUBLE);
01699 bindvar("ag",ag,DOUBLE);
01700 bindvar("b8c",b8c,DOUBLE);
01701 bindvar("d8d",d8d,DOUBLE);
01702
01703 }
01704
01705 void MEAMFrame::potential()
01706 {
01707 refreshneighborlist();
01708 screen();
01709 dscreen();
01710 rhoMEAM();
01711 kraMEAM();
01712 if(noscr==0) dscrfor();
01713 }
01714
01715 void MEAMFrame::initvars()
01716 {
01717 MDPARALLELFrame::initvars();
01718
01719 betas[0]=betas[1]=betas[2]=betas[3]=0;
01720 ts[0]=ts[1]=ts[2]=ts[3]=0;
01721 strcpy(potfile,"meamdata");
01722
01723 strcpy(lat,"unknown");
01724 strcpy(elt,"unknown");
01725 }
01726 void MEAMFrame::initparser()
01727 {
01728 MDPARALLELFrame::initparser();
01729
01730 bindvar("xzbl",&xzbl,DOUBLE);
01731 bindvar("xzbl0",&xzbl0,DOUBLE);
01732 bindvar("enable_zbl_fdimer",&enable_zbl_fdimer,INT);
01733 bindvar("enable_square_rscrn",&enable_square_rscrn,INT);
01734 }
01735
01736 int MEAMFrame::exec(char *name)
01737 {
01738 if(MDPARALLELFrame::exec(name)==0) return 0;
01739 bindcommand(name,"readMEAM",readMEAM());
01740 bindcommand(name,"printpairpot",printpairpot());
01741
01742 #ifdef _PARALLEL
01743 bindcommand(name,"Broadcast_MEAM_Param",Broadcast_MEAM_Param());
01744 #endif
01745 return -1;
01746 }
01747
01748
01749 void MEAMFrame::printpairpot()
01750 {
01751 int elti, eltj;
01752 double rmin, dr, rmax, r, phi, phip;
01753 char pairfile[200];
01754 FILE *fp;
01755
01756 elti = (int) input[0];
01757 eltj = (int) input[1];
01758 rmin = input[2];
01759 dr = input[3];
01760 rmax = input[4];
01761
01762 strcpy(pairfile,"pairpot_baskes.dat");
01763 fp=fopen(pairfile,"w");
01764 if(fp==NULL)
01765 {
01766 FATAL("printpairpot: file "<<pairfile<<" open failure.");
01767 }
01768 if(elti!=eltj)
01769 {
01770 FATAL("elti ("<<elti<<") != eltj("<<eltj<<") not implemented");
01771 }
01772 if((rmax<rmin)||(dr<=0))
01773 FATAL("rmax cannot be smaller than rmin, dr must be positive");
01774
01775 phip = 0;
01776 for(r=rmin;r<=rmax;r+=dr)
01777 {
01778 phi = phiid(r,elti);
01779 fprintf(fp,"%21.14e %21.14e %21.14e\n",r,phi,phip);
01780 }
01781 fclose(fp);
01782 INFO("results written to "<<pairfile);
01783 }
01784
01785
01786 #ifdef _PARALLEL
01787 void MEAMFrame::Broadcast_MEAM_Param()
01788 {
01789 double *buf_double;
01790 int nparam;
01791
01792
01793 if(myDomain==0) Master_to_Slave("Broadcast_MEAM_Param");
01794
01795
01796 nparam = 36;
01797 buf_double = (double *) malloc(nparam*sizeof(double));
01798 if(myDomain==0)
01799 {
01800 buf_double[0] = zsmeam;
01801 buf_double[1] = alphas;
01802 buf_double[2] = betas[0];
01803 buf_double[3] = betas[1];
01804 buf_double[4] = betas[2];
01805 buf_double[5] = betas[3];
01806 buf_double[6] = esubs;
01807 buf_double[7] = asubs;
01808 buf_double[8] = ts[0];
01809 buf_double[9] = ts[1];
01810 buf_double[10] = ts[2];
01811 buf_double[11] = ts[3];
01812 buf_double[12] = rozros;
01813 buf_double[13] = rcutmeam;
01814 buf_double[14] = cmin;
01815 buf_double[15] = cmax;
01816 buf_double[16] = repuls;
01817 buf_double[17] = attrac;
01818 buf_double[18] = legend;
01819 buf_double[19] = (double) ibarr;
01820 buf_double[20] = (double) noscr;
01821 buf_double[21] = res;
01822 buf_double[22] = alat;
01823 buf_double[23] = ielement;
01824 buf_double[24] = rcutmeam2;
01825 buf_double[25] = sconst;
01826 buf_double[26] = scnres;
01827 buf_double[27] = xzbl;
01828 buf_double[28] = xzbl0;
01829 buf_double[29] = hmeam;
01830 buf_double[30] = _RLIST;
01831 buf_double[31] = _SKIN;
01832 buf_double[32] = rcutmeam2;
01833 buf_double[33] = res;
01834 buf_double[34] = _NNM;
01835 buf_double[35] = _NIC;
01836 }
01837
01838 MPI_Bcast(buf_double,nparam,MPI_DOUBLE,0,MPI_COMM_WORLD);
01839
01840 if(myDomain!=0)
01841 {
01842
01843
01844
01845
01846
01847
01848 zsmeam = buf_double[0];
01849 alphas = buf_double[1];
01850 betas[0] = buf_double[2];
01851 betas[1] = buf_double[3];
01852 betas[2] = buf_double[4];
01853 betas[3] = buf_double[5];
01854 esubs = buf_double[6];
01855 asubs = buf_double[7];
01856 ts[0] = buf_double[8];
01857 ts[1] = buf_double[9];
01858 ts[2] = buf_double[10];
01859 ts[3] = buf_double[11];
01860 rozros = buf_double[12];
01861 rcutmeam = buf_double[13];
01862 cmin = buf_double[14];
01863 cmax = buf_double[15];
01864 repuls = buf_double[16];
01865 attrac = buf_double[17];
01866 legend = buf_double[18];
01867 ibarr = (int) buf_double[19];
01868 noscr = (int) buf_double[20];
01869 res = buf_double[21];
01870 alat = buf_double[22];
01871 ielement = buf_double[23];
01872 rcutmeam2 = buf_double[24];
01873 sconst = buf_double[25];
01874 scnres = buf_double[26];
01875 xzbl = buf_double[27];
01876 xzbl0 = buf_double[28];
01877 hmeam = buf_double[29];
01878 _RLIST = buf_double[30];
01879 _SKIN = buf_double[31];
01880 rcutmeam2 = buf_double[32];
01881 res = buf_double[33];
01882 _NNM = (int) buf_double[34];
01883 _NIC = (int) buf_double[35];
01884 }
01885
01886 free(buf_double);
01887
01888
01889 MPI_Bcast(elt,10,MPI_CHAR,0,MPI_COMM_WORLD);
01890 MPI_Bcast(lat,10,MPI_CHAR,0,MPI_COMM_WORLD);
01891 }
01892 #endif//_PARALLEL
01893
01894 #ifdef _TEST
01895
01896
01897 class MEAMFrame sim;
01898
01899
01900 #include "main.cpp"
01901
01902 #endif//_TEST
01903
01904
01905
01906
01907
01908
01909
01910
01911
01912
01913
01914
01915
01916
01917
01918
01919
01920
01921
01922
01923
01924
01925
01926
01927
01928
01929
01930
01931
01932
01933
01934
01935
01936
01937
01938
01939
01940
01941
01942
01943
01944
01945
01946
01947
01948
01949
01950
01951
01952
01953
01954
01955
01956
01957
01958
01959