void neb_sendrecv_pos(void) { int i, k, n, cpu_l, cpu_r; MPI_Status status; /* fill pos array */ for (k=0; k<NCELLS; k++) { cell *p = CELLPTR(k); for (i=0; i<p->n; i++) { n = NUMMER(p,i); pos X(n) = ORT(p,i,X); pos Y(n) = ORT(p,i,Y); pos Z(n) = ORT(p,i,Z); } } /* ranks of left/right cpus */ cpu_l = (0 == myrank) ? MPI_PROC_NULL : myrank - 1; cpu_r = (neb_nrep - 1 == myrank) ? MPI_PROC_NULL : myrank + 1; /* send positions to right, receive from left */ MPI_Sendrecv(pos, DIM*natoms, REAL, cpu_r, BUFFER_TAG, pos_l, DIM*natoms, REAL, cpu_l, BUFFER_TAG, MPI_COMM_WORLD, &status ); /* send positions to left, receive from right */ MPI_Sendrecv(pos, DIM*natoms, REAL, cpu_l, BUFFER_TAG, pos_r, DIM*natoms, REAL, cpu_r, BUFFER_TAG, MPI_COMM_WORLD, &status ); }
void init_extpot(void) { long tmpvec1[3], tmpvec2[3]; int i,k,sort; nactive_vect[0]=0; nactive_vect[1]=0; nactive_vect[2]=0; if (0==myid) { printf( "EXTPOT: choosen potential ep_key = %d\n", ep_key); printf( "EXTPOT: number of indenters = %d\n", ep_n); printf( "EXTPOT: number of walls = %d\n", ep_n-ep_nind); printf( "EXTPOT: external potential constant = %f\n",ep_a ); printf( "EXTPOT: cutoff radius = %f\n", ep_rcut ); for (i=0; i<ep_n; i++) { printf("EXTPOT: ep_pos #%d %.10g %.10g %.10g \n", i, ep_pos[i].x, ep_pos[i].y, ep_pos[i].z); printf("EXTPOT: ep_vel #%d %.10g %.10g %.10g \n", i, ep_vel[i].x, ep_vel[i].y, ep_vel[i].z); printf("EXTPOT: ep_dir #%d %.10g %.10g %.10g \n", i, ep_dir[i].x, ep_dir[i].y, ep_dir[i].z); } } /* needed if the indenter impulse should be balanced */ /* might not work with 2D, epitax, clone, superatom */ if(ep_key==1 || ep_key==2 ){ for (k=0; k<NCELLS; ++k) { cell *p = CELLPTR(k); for (i=0; i<p->n; ++i) { sort = VSORTE(p,i); nactive_vect[0] += (restrictions + sort)->x; nactive_vect[1] += (restrictions + sort)->y; nactive_vect[2] += (restrictions + sort)->z; } } #ifdef MPI tmpvec1[0] = nactive_vect[0] ; tmpvec1[1] = nactive_vect[1] ; tmpvec1[2] = nactive_vect[2] ; MPI_Allreduce( tmpvec1, tmpvec2, 3, MPI_LONG, MPI_SUM, MPI_COMM_WORLD); nactive_vect[0] = tmpvec2[0]; nactive_vect[1] = tmpvec2[1]; nactive_vect[2] = tmpvec2[2]; #endif if (0==myid) { if (nactive_vect[0] == 0 || nactive_vect[1] == 0 || nactive_vect[2] == 0) error ("ep_key=1 requires atoms free to move in all directions\n"); printf("EXTPOT: active degrees of freedom: x %ld y %ld z %ld\n",nactive_vect[0],nactive_vect[1],nactive_vect[2]); } #ifdef RELAX printf( "EXTPOT: max number of relaxation steps = %d\n", ep_max_int); printf( "EXTPOT: ekin_threshold = %f\n", glok_ekin_threshold); printf( "EXTPOT: fnorm_threshold = %f\n", fnorm_threshold); #endif } }
void add_positions(void) { int k; for (k=0; k<NCELLS; k++) { int i; cell* p; p = CELLPTR(k); for (i=0; i<p->n; i++) { AV_POS(p,i,X) += ORT(p,i,X) + SHEET(p,i,X); AV_POS(p,i,Y) += ORT(p,i,Y) + SHEET(p,i,Y); #ifndef TWOD AV_POS(p,i,Z) += ORT(p,i,Z) + SHEET(p,i,Z); #endif AV_EPOT(p,i) += POTENG(p,i); } } #ifdef NPT av_box_x.x += box_x.x; av_box_x.y += box_x.y; av_box_y.x += box_y.x; av_box_y.y += box_y.y; #ifndef TWOD av_box_x.z += box_x.z; av_box_y.z += box_y.z; av_box_z.x += box_z.x; av_box_z.y += box_z.y; av_box_z.z += box_z.z; #endif #endif avpos_cnt++; }
void pack_fcs(void) { int k, n, m, i; /* (re-)allocate data structures if necessary */ nloc = 0; for (k=0; k<NCELLS; ++k) nloc += CELLPTR(k)->n; if (nloc_max < nloc) { nloc_max = (int) (1.1 * nloc); free(pos); pos = (fcs_float*) malloc(DIM * nloc_max*sizeof(fcs_float)); free(chg); chg = (fcs_float*) malloc( nloc_max*sizeof(fcs_float)); free(field); field = (fcs_float*) malloc(DIM * nloc_max*sizeof(fcs_float)); free(pot); pot = (fcs_float*) malloc( nloc_max*sizeof(fcs_float)); if ((NULL==pos) || (NULL==chg) || (NULL==field) || (NULL==pot)) error("Could not allocate fcs data"); } /* if we do only fcs, we have to clear forces and energies */ #ifdef PAIR if (!have_potfile && !have_pre_pot) #endif clear_forces(); /* collect data from cell array */ n=0; m=0; for (k=0; k<NCELLS; ++k) { cell *p = CELLPTR(k); for (i=0; i<p->n; ++i) { pos[n++] = ORT(p,i,X); pos[n++] = ORT(p,i,Y); pos[n++] = ORT(p,i,Z); chg[m++] = CHARGE(p,i); } } /* apply periodic boundaries */ do_boundaries_fcs(nloc, pos); }
void add_presstensors(void) { int k; for (k=0; k<NCELLS; k++) { int i; cell* p; p = CELLPTR(k); for (i=0; i<p->n; i++) { AVPRESSTENS(p,i,xx) += PRESSTENS(p,i,xx); AVPRESSTENS(p,i,yy) += PRESSTENS(p,i,yy); AVPRESSTENS(p,i,xy) += PRESSTENS(p,i,xy); #ifndef TWOD AVPRESSTENS(p,i,zz) += PRESSTENS(p,i,zz); AVPRESSTENS(p,i,zx) += PRESSTENS(p,i,zx); AVPRESSTENS(p,i,yz) += PRESSTENS(p,i,yz); #endif } } }
void maxwell(real temp) { int k; vektor tot_impuls; int nactive_x, nactive_y, nactive_z; static long dummy = 0; int slice; #ifdef DAMP real tmp1,tmp2,tmp3,f,maxax,maxax2; #endif #ifdef LASER real depth; #endif #ifdef UNIAX real xisq; real xi0; real xi1; real xi2; real dot; real norm; real osq; #endif real TEMP; real scale, xx, tmp; int num, nhalf, typ; TEMP = temp; tot_impuls.x = 0.0; nactive_x = 0; tot_impuls.y = 0.0; nactive_y = 0; #ifndef TWOD tot_impuls.z = 0.0; nactive_z = 0; #endif /* set temperature */ for (k=0; k<NCELLS; ++k) { int i; cell *p; vektor *rest; p = CELLPTR(k); for (i=0; i<p->n; ++i) { #ifdef LASER depth = laser_dir.x * ORT(p,i,X) + laser_dir.y * ORT(p,i,Y) #ifndef TWOD + laser_dir.z * ORT(p,i,Z) #endif ; depth -= laser_offset; if (depth < 0) { depth=0; /* we don't want to exceed laser_delta_temp */ } TEMP = laser_delta_temp * exp(-laser_mu * depth); TEMP += temperature; /* add base temperature of sample */ #endif /* LASER */ #ifdef DAMP /* Calculate stadium function f */ maxax = MAX(MAX(stadium.x,stadium.y),stadium.z); maxax2 = MAX(MAX(stadium2.x,stadium2.y),stadium2.z); tmp1 = (stadium2.x == 0) ? 0 : SQR((ORT(p,i,X)-center.x)/(2.0*stadium2.x)); tmp2 = (stadium2.y == 0) ? 0 : SQR((ORT(p,i,Y)-center.y)/(2.0*stadium2.y)); tmp3 = (stadium2.z == 0) ? 0 : SQR((ORT(p,i,Z)-center.z)/(2.0*stadium2.z)); f = (tmp1+tmp2+tmp3-SQR(maxax/(2.0*maxax2)))/\ (.25- SQR(maxax/(2.0*maxax2))); if (f<= 0.0) f = 0.0; else if (f>1.0) f = 1.0; /* we smooth the stadium function: to get a real bath tub ! fully damped atoms have temp=0, temperature gradient determined by bath tub */ if (f !=0) TEMP = damptemp * (1.0 - 0.5 * (1 + sin(-M_PI/2.0 + M_PI*f))); else TEMP= temp; #endif #ifdef FTG /* calc slice and set TEMP */ tmp = ORT(p,i,X) / box_x.x; slice = (int) nslices * tmp; if (slice<0) slice = 0; if (slice>=nslices) slice = nslices -1;; TEMP= Tleft + (Tright-Tleft)*(slice-nslices_Left+1) / (real) (nslices-nslices_Left-nslices_Right+1); if (slice>=nslices-nslices_Right) TEMP = Tright; if (slice<nslices_Left) TEMP= Tleft; #endif tmp = sqrt(TEMP * MASSE(p,i)); rest = restrictions + VSORTE(p,i); #ifndef RIGID IMPULS(p,i,X) = imd_gaussian(tmp) * rest->x; IMPULS(p,i,Y) = imd_gaussian(tmp) * rest->y; #ifndef TWOD IMPULS(p,i,Z) = imd_gaussian(tmp) * rest->z; #endif #else /* superatoms get velocity zero */ if (superatom[VSORTE(p,i)]>-1) { IMPULS(p,i,X) = 0.0; IMPULS(p,i,Y) = 0.0; #ifndef TWOD IMPULS(p,i,Z) = 0.0; #endif } #endif nactive_x += (int) rest->x; nactive_y += (int) rest->y; #ifndef TWOD nactive_z += (int) rest->z; #endif tot_impuls.x += IMPULS(p,i,X); tot_impuls.y += IMPULS(p,i,Y); #ifndef TWOD tot_impuls.z += IMPULS(p,i,Z); #endif #ifdef UNIAX /* angular velocities for uniaxial molecules */ /* choose a random vector in space */ do { xi1 = 2.0 * drand48() - 1.0 ; xi2 = 2.0 * drand48() - 1.0 ; xisq = xi1 * xi1 + xi2 * xi2 ; } while ( xisq >= 1.0 ); xi0 = sqrt( 1.0 - xisq ) ; DREH_IMPULS(p,i,X) = 2.0 * xi1 * xi0 ; DREH_IMPULS(p,i,Y) = 2.0 * xi2 * xi0 ; DREH_IMPULS(p,i,Z) = 1.0 - 2.0 * xisq ; /* constrain the vector to be perpendicular to the molecule */ dot = SPRODN(DREH_IMPULS,p,i,ACHSE,p,i); DREH_IMPULS(p,i,X) -= dot * ACHSE(p,i,X) ; DREH_IMPULS(p,i,Y) -= dot * ACHSE(p,i,Y) ; DREH_IMPULS(p,i,Z) -= dot * ACHSE(p,i,Z) ; /* renormalize vector */ osq = SPRODN(DREH_IMPULS,p,i,DREH_IMPULS,p,i); norm = SQRT( osq ); DREH_IMPULS(p,i,X) /= norm ; DREH_IMPULS(p,i,Y) /= norm ; DREH_IMPULS(p,i,Z) /= norm ; /* choose the magnitude of the angular momentum */ osq = - 2.0 * uniax_inert * TEMP * log( drand48() ) ; norm = sqrt( osq ); DREH_IMPULS(p,i,X) *= norm ; DREH_IMPULS(p,i,Y) *= norm ; DREH_IMPULS(p,i,Z) *= norm ; #endif /* UNIAX */ #ifdef SHOCK /* plate against bulk */ if (shock_mode == 1) { if ( ORT(p,i,X) < shock_strip ) IMPULS(p,i,X) += shock_speed * MASSE(p,i); } /* two halves against one another */ if (shock_mode == 2) { if ( ORT(p,i,X) < box_x.x*0.5 ) IMPULS(p,i,X) += shock_speed * MASSE(p,i); else IMPULS(p,i,X) -= shock_speed * MASSE(p,i); } /* bulk against wall */ if (shock_mode == 3) IMPULS(p,i,X) += shock_speed * MASSE(p,i); #endif } } #ifdef CLONE /* compute the total momentum afresh */ tot_impuls.x = 0.0; tot_impuls.y = 0.0; #ifndef TWOD tot_impuls.z = 0.0; #endif /* set velocities of clones equal */ for (k=0; k<NCELLS; k++) { int i, j; cell *p; p = CELLPTR(k); for (i=0; i<p->n; i+=nclones) { tot_impuls.x += nclones * IMPULS(p,i,X); tot_impuls.y += nclones * IMPULS(p,i,Y); #ifndef TWOD tot_impuls.z += nclones * IMPULS(p,i,Z); #endif for (j=1; j<nclones; j++) { IMPULS(p,i+j,X) = IMPULS(p,i,X); IMPULS(p,i+j,Y) = IMPULS(p,i,Y); #ifndef TWOD IMPULS(p,i+j,Z) = IMPULS(p,i,Z); #endif } } } #endif /* CLONE */ tot_impuls.x = nactive_x == 0 ? 0.0 : tot_impuls.x / nactive_x; tot_impuls.y = nactive_y == 0 ? 0.0 : tot_impuls.y / nactive_y; #ifndef TWOD tot_impuls.z = nactive_z == 0 ? 0.0 : tot_impuls.z / nactive_z; #endif /* correct center of mass momentum */ for (k=0; k<NCELLS; ++k) { int i; cell *p; vektor *rest; p = CELLPTR(k); for (i=0; i<p->n; ++i) { rest = restrictions + VSORTE(p,i); IMPULS(p,i,X) -= tot_impuls.x * rest->x; IMPULS(p,i,Y) -= tot_impuls.y * rest->y; #ifndef TWOD IMPULS(p,i,Z) -= tot_impuls.z * rest->z; #endif } } }
void calc_extpot(void) { int k, i, n; int isinx,isiny,isinz; real tmpvec1[4], tmpvec2[4]; vektor d,addforce,totaddforce; real dd,cc; real dn,ddn,ee; vektor force; real tmp_virial; #ifdef P_AXIAL vektor tmp_vir_vect; #endif real pot_zwi, pot_grad; int col, is_short=0; for (k=0; k<ep_n; k++) { ep_fext[k] = 0.0; ep_xmax[k] = 0.0; ep_ymax[k] = 0.0; ep_atomsincontact[k]=0; ep_xmin[k] = 1.e8; ep_ymin[k] = 1.e8; } if(ep_key == 0) { /* default: original harmonic potential */ for (k=0; k<NCELLS; ++k) { cell *p = CELLPTR(k); for (i=0; i<p->n; ++i) { for (n=0; n<ep_n; ++n) { isinx= ep_dir[n].x; isiny= ep_dir[n].y; isinz= ep_dir[n].z; d.x = ep_pos[n].x - ORT(p,i,X); d.y = ep_pos[n].y - ORT(p,i,Y); d.z = ep_pos[n].z - ORT(p,i,Z); dn = SPROD(d,ep_dir[n]); /* spherical indentor*/ if (n<ep_nind) { if (dn > -ep_rcut) { real d2 = SPROD(d,d); real d1 = SQRT(d2); dd = ep_rcut - d1; if (dd > 0.0) { real f = ep_a * dd * dd / d1; /* force on atoms and indentor */ KRAFT(p,i,X) -= f * d.x; KRAFT(p,i,Y) -= f * d.y; KRAFT(p,i,Z) -= f * d.z; ep_fext[n] += f * ABS(dn); /* normal force on indentor */ ep_atomsincontact[n]++; /* for determination of contact area */ if(isinz) { ep_xmax[n] = MAX(ep_xmax[n], ORT(p,i,X) ); ep_ymax[n] = MAX(ep_ymax[n], ORT(p,i,Y) ); ep_xmin[n] = MIN(ep_xmin[n], ORT(p,i,X) ); ep_ymin[n] = MIN(ep_ymin[n], ORT(p,i,Y) ); } else if(isiny) { ep_xmax[n] = MAX(ep_xmax[n], ORT(p,i,X) ); ep_ymax[n] = MAX(ep_ymax[n], ORT(p,i,Z) ); ep_xmin[n] = MIN(ep_xmin[n], ORT(p,i,X) ); ep_ymin[n] = MIN(ep_ymin[n], ORT(p,i,Z) ); } else { ep_xmax[n] = MAX(ep_xmax[n], ORT(p,i,Y) ); ep_ymax[n] = MAX(ep_ymax[n], ORT(p,i,Z) ); ep_xmin[n] = MIN(ep_xmin[n], ORT(p,i,Y) ); ep_ymin[n] = MIN(ep_ymin[n], ORT(p,i,Z) ); } } } } /* potential wall */ else { if (dn*dn < ep_rcut*ep_rcut) { real d1 = (dn>0) ? dn : -1*dn ; dd = ep_rcut - d1; if (dd > 0.0) { ep_atomsincontact[n]++; real f = ep_a * dd * dd / d1; /* force on atoms and indentor */ KRAFT(p,i,X) += f * ep_dir[n].x; KRAFT(p,i,Y) += f * ep_dir[n].y; KRAFT(p,i,Z) += f * ep_dir[n].z; ep_fext[n] += f; /* magnitude of force on wall */ } } } } } } } else if(ep_key == 1) /* Ju Li's spherical indenter, see PRB 67, 104105 */ { totaddforce.x=0.0; totaddforce.y=0.0; totaddforce.z=0.0; for (k=0; k<NCELLS; ++k) { cell *p = CELLPTR(k); for (i=0; i<p->n; ++i) { for (n=0; n<ep_n; ++n) { isinx= ep_dir[n].x; isiny= ep_dir[n].y; isinz= ep_dir[n].z; if (n<ep_nind) { d.x = ORT(p,i,X)-ep_pos[n].x; d.y = ORT(p,i,Y)-ep_pos[n].y; d.z = ORT(p,i,Z)-ep_pos[n].z; dn = SPROD(d,ep_dir[n]); dd = SPROD(d,d); if ( dd < ep_rcut*ep_rcut) { ep_atomsincontact[n]++; /* for the determination of the contact area */ if(isinz) { ep_xmax[n] = MAX(ep_xmax[n], ORT(p,i,X) ); ep_ymax[n] = MAX(ep_ymax[n], ORT(p,i,Y) ); ep_xmin[n] = MIN(ep_xmin[n], ORT(p,i,X) ); ep_ymin[n] = MIN(ep_ymin[n], ORT(p,i,Y) ); } else if(isiny) { ep_xmax[n] = MAX(ep_xmax[n], ORT(p,i,X) ); ep_ymax[n] = MAX(ep_ymax[n], ORT(p,i,Z) ); ep_xmin[n] = MIN(ep_xmin[n], ORT(p,i,X) ); ep_ymin[n] = MIN(ep_ymin[n], ORT(p,i,Z) ); } else { ep_xmax[n] = MAX(ep_xmax[n], ORT(p,i,Y) ); ep_ymax[n] = MAX(ep_ymax[n], ORT(p,i,Z) ); ep_xmin[n] = MIN(ep_xmin[n], ORT(p,i,Y) ); ep_ymin[n] = MIN(ep_ymin[n], ORT(p,i,Z) ); } if(have_extpotfile == 1){ PAIR_INT3(pot_zwi, pot_grad, ext_pot, n, ep_nind, dd, is_short); tot_pot_energy += pot_zwi; force.x = -1.0* pot_grad * d.x; force.y = -1.0* pot_grad * d.y; force.z = -1.0* pot_grad * d.z; KRAFT(p,i,X) += force.x; KRAFT(p,i,Y) += force.y; KRAFT(p,i,Z) += force.z; totaddforce.x += force.x; totaddforce.y += force.y; totaddforce.z += force.z; ep_fext[n] += -pot_grad * ABS(dn); /* normal force on indentor */ #ifdef P_AXIAL tmp_vir_vect.x -= d.x * force.x; tmp_vir_vect.y -= d.y * force.y; #ifndef TWOD tmp_vir_vect.z -= d.z * force.z; #endif #else tmp_virial -= dd * pot_grad; #endif #ifdef STRESS_TENS if (do_press_calc) { PRESSTENS(p,i,xx) -= d.x * force.x; PRESSTENS(p,i,yy) -= d.y * force.y; PRESSTENS(p,i,xy) -= d.x * force.y; #ifndef TWOD PRESSTENS(p,i,zz) -= d.z * force.z; PRESSTENS(p,i,yz) -= d.y * force.z; PRESSTENS(p,i,zx) -= d.z * force.x; #endif } #endif } /* old version of extpot, kept for downwards compatibility */ else{ ddn= sqrt(dd); cc = (ep_rcut - ddn)/ep_a; if (cc > UPPER_EXP) cc = UPPER_EXP; if (cc < LOWER_EXP) cc = LOWER_EXP; ee = exp(cc - 1.0/cc); tot_pot_energy += ee; POTENG(p,i) += ee; ee = ee / ep_a / ddn * (1.0 + 1.0 /(cc*cc)); KRAFT(p,i,X) += ee * d.x; KRAFT(p,i,Y) += ee * d.y; KRAFT(p,i,Z) += ee * d.z; totaddforce.x += ee * d.x; totaddforce.y += ee * d.y; totaddforce.z += ee * d.z; ep_fext[n] += ee * ABS(dn); /* normal force on indentor */ } } } } } } #ifdef MPI tmpvec1[0] = totaddforce.x ; tmpvec1[1] = totaddforce.y ; tmpvec1[2] = totaddforce.z ; // printf("before totaddforcereduce allreduce\n");fflush(stdout); MPI_Allreduce( tmpvec1, tmpvec2, 4, REAL, MPI_SUM, cpugrid); // printf("after totaddforce allreduce\n");fflush(stdout); totaddforce.x = tmpvec2[0]; totaddforce.y = tmpvec2[1]; totaddforce.z = tmpvec2[2]; #endif /* no need for a wall as the total additional impuls is substracted */ totaddforce.x *= 1.0/nactive_vect[0]; totaddforce.y *= 1.0/nactive_vect[1]; totaddforce.z *= 1.0/nactive_vect[2]; for (k=0; k<NCELLS; ++k) { cell *p = CELLPTR(k); for (i=0; i<p->n; ++i) { KRAFT(p,i,X) -= totaddforce.x; KRAFT(p,i,Y) -= totaddforce.y; KRAFT(p,i,Z) -= totaddforce.z; } } } else if(ep_key == 2) /* Ju Li's spherical indenter made flat, see PRB 67, 104105 with subtraction of total additional impulse works only with indentation directions parallel to box vectors*/ { // vektor d,addforce,totaddforce; //real dd,cc; //real dn,ddn,ee; totaddforce.x=0.0; totaddforce.y=0.0; totaddforce.z=0.0; for (k=0; k<NCELLS; ++k) { cell *p = CELLPTR(k); for (i=0; i<p->n; ++i) { for (n=0; n<ep_n; ++n) { isinx= ep_dir[n].x; isiny= ep_dir[n].y; isinz= ep_dir[n].z; // vektor d; // real dn; d.x = (ep_dir[n].x==0) ? 0 : (ORT(p,i,X)-ep_pos[n].x); d.y = (ep_dir[n].y==0) ? 0 : (ORT(p,i,Y)-ep_pos[n].y); d.z = (ep_dir[n].z==0) ? 0 : (ORT(p,i,Z)-ep_pos[n].z); dn = SPROD(d,ep_dir[n]); dd = SPROD(d,d); if ( dd < ep_rcut*ep_rcut) { /* for the determination of the contact area */ ep_atomsincontact[n]++; if(isinz) { ep_xmax[n] = MAX(ep_xmax[n], ORT(p,i,X) ); ep_ymax[n] = MAX(ep_ymax[n], ORT(p,i,Y) ); ep_xmin[n] = MIN(ep_xmin[n], ORT(p,i,X) ); ep_ymin[n] = MIN(ep_ymin[n], ORT(p,i,Y) ); } else if(isiny) { ep_xmax[n] = MAX(ep_xmax[n], ORT(p,i,X) ); ep_ymax[n] = MAX(ep_ymax[n], ORT(p,i,Z) ); ep_xmin[n] = MIN(ep_xmin[n], ORT(p,i,X) ); ep_ymin[n] = MIN(ep_ymin[n], ORT(p,i,Z) ); } else { ep_xmax[n] = MAX(ep_xmax[n], ORT(p,i,Y) ); ep_ymax[n] = MAX(ep_ymax[n], ORT(p,i,Z) ); ep_xmin[n] = MIN(ep_xmin[n], ORT(p,i,Y) ); ep_ymin[n] = MIN(ep_ymin[n], ORT(p,i,Z) ); } if(have_extpotfile == 1){ PAIR_INT3(pot_zwi, pot_grad, ext_pot, n, ep_nind, dd, is_short); tot_pot_energy += pot_zwi; force.x = -1.0* pot_grad * d.x; force.y = -1.0* pot_grad * d.y; force.z = -1.0* pot_grad * d.z; KRAFT(p,i,X) += force.x; KRAFT(p,i,Y) += force.y; KRAFT(p,i,Z) += force.z; totaddforce.x += force.x; totaddforce.y += force.y; totaddforce.z += force.z; ep_fext[n] += -pot_grad * ABS(dn); /* normal force on indentor */ #ifdef P_AXIAL tmp_vir_vect.x -= d.x * force.x; tmp_vir_vect.y -= d.y * force.y; #ifndef TWOD tmp_vir_vect.z -= d.z * force.z; #endif #else tmp_virial -= dd * pot_grad; #endif #ifdef STRESS_TENS if (do_press_calc) { PRESSTENS(p,i,xx) -= d.x * force.x; PRESSTENS(p,i,yy) -= d.y * force.y; PRESSTENS(p,i,xy) -= d.x * force.y; #ifndef TWOD PRESSTENS(p,i,zz) -= d.z * force.z; PRESSTENS(p,i,yz) -= d.y * force.z; PRESSTENS(p,i,zx) -= d.z * force.x; #endif } #endif } /* old version of extpot, kept for downwards compatibility */ else{ ddn= sqrt(dd); cc = (ep_rcut - ddn)/ep_a; if (cc > UPPER_EXP) cc = UPPER_EXP; if (cc < LOWER_EXP) cc = LOWER_EXP; ee = exp(cc - 1.0/cc); tot_pot_energy += ee; POTENG(p,i) += ee; ee = ee / ep_a / ddn * (1.0 + 1.0 /(cc*cc)); KRAFT(p,i,X) += ee * d.x; KRAFT(p,i,Y) += ee * d.y; KRAFT(p,i,Z) += ee * d.z; totaddforce.x += ee * d.x; totaddforce.y += ee * d.y; totaddforce.z += ee * d.z; ep_fext[n] += ee * ABS(dn); /* normal force on indentor */ } } } } } #ifdef MPI tmpvec1[0] = totaddforce.x ; tmpvec1[1] = totaddforce.y ; tmpvec1[2] = totaddforce.z ; // printf("before totaddforcereduce allreduce\n");fflush(stdout); MPI_Allreduce( tmpvec1, tmpvec2, 4, REAL, MPI_SUM, cpugrid); // printf("after totaddforce allreduce\n");fflush(stdout); totaddforce.x = tmpvec2[0]; totaddforce.y = tmpvec2[1]; totaddforce.z = tmpvec2[2]; #endif /* no need for a wall as the total additional impuls is substracted */ totaddforce.x *= 1.0/nactive_vect[0]; totaddforce.y *= 1.0/nactive_vect[1]; totaddforce.z *= 1.0/nactive_vect[2]; for (k=0; k<NCELLS; ++k) { cell *p = CELLPTR(k); for (i=0; i<p->n; ++i) { KRAFT(p,i,X) -= totaddforce.x; KRAFT(p,i,Y) -= totaddforce.y; KRAFT(p,i,Z) -= totaddforce.z; } } } else if(ep_key == 3) /* Ju Li's spherical indenter made flat, see PRB 67, 104105 without subtraction of the total additional impulse works only with indentation directions parallel to box vectors*/ { // vektor d,addforce,totaddforce; // real dd,cc; // real dn,ddn,ee; totaddforce.x=0.0; totaddforce.y=0.0; totaddforce.z=0.0; for (k=0; k<NCELLS; ++k) { cell *p = CELLPTR(k); for (i=0; i<p->n; ++i) { for (n=0; n<ep_n; ++n) { isinx= ep_dir[n].x; isiny= ep_dir[n].y; isinz= ep_dir[n].z; vektor d; real dn; d.x = (ep_dir[n].x==0) ? 0 : (ORT(p,i,X)-ep_pos[n].x) ; d.y = (ep_dir[n].y==0) ? 0 : (ORT(p,i,Y)-ep_pos[n].y); d.z = (ep_dir[n].z==0) ? 0 : (ORT(p,i,Z)-ep_pos[n].z); dn = SPROD(d,ep_dir[n]); dd = SPROD(d,d); if ( dd < ep_rcut*ep_rcut) { /* for the determination of the contact area */ ep_atomsincontact[n]++; if(isinz) { ep_xmax[n] = MAX(ep_xmax[n], ORT(p,i,X) ); ep_ymax[n] = MAX(ep_ymax[n], ORT(p,i,Y) ); ep_xmin[n] = MIN(ep_xmin[n], ORT(p,i,X) ); ep_ymin[n] = MIN(ep_ymin[n], ORT(p,i,Y) ); } else if(isiny) { ep_xmax[n] = MAX(ep_xmax[n], ORT(p,i,X) ); ep_ymax[n] = MAX(ep_ymax[n], ORT(p,i,Z) ); ep_xmin[n] = MIN(ep_xmin[n], ORT(p,i,X) ); ep_ymin[n] = MIN(ep_ymin[n], ORT(p,i,Z) ); } else { ep_xmax[n] = MAX(ep_xmax[n], ORT(p,i,Y) ); ep_ymax[n] = MAX(ep_ymax[n], ORT(p,i,Z) ); ep_xmin[n] = MIN(ep_xmin[n], ORT(p,i,Y) ); ep_ymin[n] = MIN(ep_ymin[n], ORT(p,i,Z) ); } if(have_extpotfile == 1){ PAIR_INT3(pot_zwi, pot_grad, ext_pot, n, ep_nind, dd, is_short); tot_pot_energy += pot_zwi; force.x = - pot_grad * d.x; force.y = - pot_grad * d.y; force.z = - pot_grad * d.z; KRAFT(p,i,X) += force.x; KRAFT(p,i,Y) += force.y; KRAFT(p,i,Z) += force.z; totaddforce.x += force.x; totaddforce.y += force.y; totaddforce.z += force.z; ep_fext[n] += -pot_grad * ABS(dn); /* normal force on indentor */ #ifdef P_AXIAL tmp_vir_vect.x -= d.x * force.x; tmp_vir_vect.y -= d.y * force.y; #ifndef TWOD tmp_vir_vect.z -= d.z * force.z; #endif #else tmp_virial -= dd * pot_grad; #endif #ifdef STRESS_TENS if (do_press_calc) { PRESSTENS(p,i,xx) -= d.x * force.x; PRESSTENS(p,i,yy) -= d.y * force.y; PRESSTENS(p,i,xy) -= d.x * force.y; #ifndef TWOD PRESSTENS(p,i,zz) -= d.z * force.z; PRESSTENS(p,i,yz) -= d.y * force.z; PRESSTENS(p,i,zx) -= d.z * force.x; #endif } #endif } /* old version of extpot, kept for downwards compatibility */ else{ ddn= sqrt(dd); cc = (ep_rcut - ddn)/ep_a; if (cc > UPPER_EXP) cc = UPPER_EXP; if (cc < LOWER_EXP) cc = LOWER_EXP; ee = exp(cc - 1.0/cc); tot_pot_energy += ee; POTENG(p,i) += ee; ee = ee / ep_a / ddn * (1.0 + 1.0 /(cc*cc)); KRAFT(p,i,X) += ee * d.x; KRAFT(p,i,Y) += ee * d.y; KRAFT(p,i,Z) += ee * d.z; ep_fext[n] += ee * ABS(dn); /* normal force on indentor */ } } } } } } else { error("Error: external potential ep_key not defined.\n"); } #ifdef P_AXIAL vir_xx += tmp_vir_vect.x; vir_yy += tmp_vir_vect.y; virial += tmp_vir_vect.x; virial += tmp_vir_vect.y; #ifndef TWOD vir_zz += tmp_vir_vect.z; virial += tmp_vir_vect.z; #endif #else virial += tmp_virial; #endif }
/* This is laser_rescale_mode == 4 */ void laser_rescale_ttm() { /* This function just writes the exponential source term into the FD net. * Heating of the electrons occurs later, * in the numerical solution of the pdeq. */ int i,j,k; real exp_gauss_time_etc, gauss_time_squared, gauss_time_squared1, depth; gauss_time_squared = timestep * steps - laser_t_0; gauss_time_squared *= gauss_time_squared; exp_gauss_time_etc = exp(-gauss_time_squared/laser_sigma_t_squared/2.0) * laser_p_peak ;/* * fd_h.x*fd_h.y*fd_h.z not needed */ if (laser_t_1>0) { gauss_time_squared1 = timestep * steps - laser_t_1; gauss_time_squared1 *= gauss_time_squared1; exp_gauss_time_etc += exp(-gauss_time_squared1/laser_sigma_t1_squared/2.0) * laser_p_peak1;/* * fd_h.x*fd_h.y*fd_h.z not needed */ } /* loop over all FD cells, excluding ghost layers */ for(i=1; i<local_fd_dim.x-1; ++i) { for(j=1; j<local_fd_dim.y-1; ++j) { for(k=1; k<local_fd_dim.z-1; ++k) { depth = ttm_calc_depth(i,j,k); l2[i][j][k].source = l1[i][j][k].source = exp(-laser_mu*depth) * exp_gauss_time_etc; } } } #ifdef PDECAY for (k=0; k<NCELLS; k++) { cell *p; p = CELLPTR(k); int i; for (i=0; i < p->n; i++) { if( ORT(p,i,X) > ramp_start ) { switch ( pdecay_mode){ case 0: { double m = 1 / ( (ramp_end - ramp_start ) ) ; double b = - ramp_start / ( ramp_end - ramp_start ); double xia = ORT(p,i,X) *m + b ; IMPULS(p,i,X) *= 1.0- ( ORT(p,i,X) *m + b) ; if(steps==-1) printf(" %f %f \n", ORT(p,i,X),1.0 - ( ORT(p,i,X) * m + b) ); break; } case 1: { double a = 1.0 / ( ramp_start - ramp_end); a *= a; IMPULS(p,i,X) *= a * ( ORT(p,i,X) - ramp_end ) * ( ORT(p,i,X) - ramp_end ); if(steps==-1) printf("%f %f \n", ORT(p,i,X), a * ( ORT(p,i,X) - ramp_end ) * ( ORT(p,i,X) - ramp_end )); break; } break; case 2: { double m = 1 / ( (ramp_end - ramp_start ) ) ; double b = - ramp_start / ( ramp_end - ramp_start ); KRAFT(p,i,X) -= ( IMPULS(p,i,X)/MASSE(p,i)) * ( ORT(p,i,X) *m + b ) * xipdecay; if(steps==-1) printf(" %f %f \n", ORT(p,i,X),ORT(p,i,X) *m + b); break; } break; case 3: { double a = 1.0 / ( ramp_end - ramp_start); a *= a; KRAFT(p,i,X) -= ( IMPULS(p,i,X)/MASSE(p,i)) * xipdecay * a * ( ORT(p,i,X) - ramp_start ) * ( ORT(p,i,X) - ramp_start ); if(steps==-1) printf("%f %f \n", ORT(p,i,X), a * ( ORT(p,i,X) - ramp_start ) * ( ORT(p,i,X) - ramp_start ) ); break; } break; default: error("Illegal value for parameter pdecay_mode.\n"); break; } } } } #endif }
void laser_rescale_2() /* Instead of just rescaling, add the momentum increment in random direction. Only thereafter rescale so the right amount of energy gets absorbed. Involves several square roots, probably slower */ { real exp_gauss_time_etc, gauss_time_squared, gauss_time_squared1; int k; gauss_time_squared = timestep * steps - laser_t_0; gauss_time_squared *= gauss_time_squared; exp_gauss_time_etc = exp(-gauss_time_squared/laser_sigma_t_squared/2.0) * laser_p_peak * timestep * laser_atom_vol; if (laser_t_1>0) { gauss_time_squared1 = timestep * steps - laser_t_1; gauss_time_squared1 *= gauss_time_squared1; exp_gauss_time_etc += exp(-gauss_time_squared1/laser_sigma_t1_squared/2.0) * laser_p_peak1 * timestep * laser_atom_vol; } for (k=0; k<NCELLS; k++) { cell *p; p = CELLPTR(k); int i; for (i=0; i < p->n; i++) { real p_0_square; /* square of initial atom momentum */ real p_0; real de; /* Kinetic energy to be added to the atom */ real dp; real depth; /* Distance from origin along laser_dir */ real tmpx, tmpy, tmpsqr; #ifndef TWOD real tmpz; #endif real scale_p; /* Scale factor for atom momentum */ p_0_square = SPRODN(IMPULS,p,i,IMPULS,p,i); p_0=sqrt(p_0_square); depth = laser_calc_depth(p,i); #ifdef LASERYZ /* spacial dependence of laser beam */ double x = ORT(p,i,X); double y = ORT(p,i,Y); double z = ORT(p,i,Z); de = exp( -laser_mu*depth ) * exp_gauss_time_etc * laser_intensity_profile(x,y,z); #else de = exp(-laser_mu*depth) * exp_gauss_time_etc; #endif dp = sqrt( p_0_square + 2*de*MASSE(p,i) ) - p_0; /* Momentum increment is to point in a random direction... */ #ifndef TWOD rand_uvec_3d(&tmpx, &tmpy, &tmpz); #else rand_uvec_2d(&tmpx, &tmpy); #endif IMPULS(p,i,X) += tmpx * dp; IMPULS(p,i,Y) += tmpy * dp; #ifndef TWOD IMPULS(p,i,Z) += tmpz * dp; #endif /* rescale present momentum to an absolute value of dp + p_0 */ scale_p = (p_0 + dp) / SQRT( SPRODN(IMPULS,p,i,IMPULS,p,i) ); IMPULS(p,i,X) *= scale_p; IMPULS(p,i,Y) *= scale_p; #ifndef TWOD IMPULS(p,i,Z) *= scale_p; #endif } } } /* void laser_rescale_2()*/
void laser_rescale_1() { /* Rescale all atom velocities so that exactly the same amount of kinetic energy is added to the atoms at a given depth. The direction of the momentum vectors remains unchanged, except for zero-velocity atoms, which will be given a random direction.*/ real exp_gauss_time_etc, gauss_time_squared, gauss_time_squared1; int k; real cosph, sinph, theta, costh, sinth; gauss_time_squared = timestep * steps - laser_t_0; gauss_time_squared *= gauss_time_squared; exp_gauss_time_etc = exp(-gauss_time_squared/laser_sigma_t_squared/2.0) * laser_p_peak * timestep * laser_atom_vol; if (laser_t_1>0) { gauss_time_squared1 = timestep * steps - laser_t_1; gauss_time_squared1 *= gauss_time_squared1; exp_gauss_time_etc += exp(-gauss_time_squared1/laser_sigma_t1_squared/2.0) * laser_p_peak1 * timestep * laser_atom_vol; } for (k=0; k<NCELLS; k++) { cell *p; p = CELLPTR(k); int i; for (i=0; i < p->n; i++) { real p_0_square; /* square of initial atom momentum */ real de; /* Kinetic energy to be added to the atom */ real depth; /* Distance from origin along laser_dir */ real tmpx, tmpy, tmpsqr; #ifndef TWOD real tmpz; #endif real scale_p; /* Scale factor for atom momentum */ p_0_square = SPRODN(IMPULS,p,i,IMPULS,p,i); depth = laser_calc_depth(p,i); #ifdef PDECAY if( ORT(p,i,X) > ramp_start ) { switch ( pdecay_mode){ case 0: { /* y= m * x + b */ double m = 1 / ( (ramp_end - ramp_start ) ) ; double b = - ramp_start / ( ramp_end - ramp_start ); double xia = ORT(p,i,X) *m + b ; IMPULS(p,i,X) *= 1.0 - ( ORT(p,i,X) *m + b) ; if(steps==-1) printf(" %f %f \n", ORT(p,i,X),1.0 - ( ORT(p,i,X) * m + b) ); break; } case 1: { double a = 1.0 / ( ramp_start - ramp_end); a *= a; IMPULS(p,i,X) *= a * ( ORT(p,i,X) - ramp_end ) * ( ORT(p,i,X) - ramp_end ); if(steps==-1) printf("%f %f \n", ORT(p,i,X), a * ( ORT(p,i,X) - ramp_end ) * ( ORT(p,i,X) - ramp_end )); break; } break; case 2: { double m = 1 / ( (ramp_end - ramp_start ) ) ; double b = - ramp_start / ( ramp_end - ramp_start ); KRAFT(p,i,X) -= ( IMPULS(p,i,X)/MASSE(p,i)) * ( ORT(p,i,X) *m + b ) * xipdecay; if(steps==-1) printf(" %f %f \n", ORT(p,i,X),ORT(p,i,X) *m + b); break; } break; case 3: { double a = 1.0 / ( ramp_end - ramp_start); a *= a; KRAFT(p,i,X) -= ( IMPULS(p,i,X)/MASSE(p,i)) * xipdecay * a * ( ORT(p,i,X) - ramp_start ) * ( ORT(p,i,X) - ramp_start ); if(steps==-1) printf("%f %f \n", ORT(p,i,X), a * ( ORT(p,i,X) - ramp_start ) * ( ORT(p,i,X) - ramp_start ) ); break; } break; default: error("Illegal value for parameter pdecay_mode.\n"); break; } } #endif #ifdef LASERYZ /* spacial dependence of laser beam */ double x = ORT(p,i,X); double y = ORT(p,i,Y); double z = ORT(p,i,Z); de = exp( -laser_mu*depth ) * exp_gauss_time_etc * laser_intensity_profile(x,y,z); #else de = exp(-laser_mu*depth) * exp_gauss_time_etc; #endif if ( p_0_square == 0.0 ) { /* we need a direction for the momentum. */ #ifndef TWOD /* find random 3d unit vector */ rand_uvec_3d(&tmpx, &tmpy, &tmpz); #else /* find random 2d unit vector */ rand_uvec_2d(&tmpx, &tmpy); #endif scale_p = sqrt( de * 2.0 * MASSE(p,i) ); IMPULS(p,i,X) = tmpx * scale_p; IMPULS(p,i,Y) = tmpy * scale_p; #ifndef TWOD IMPULS(p,i,Z) = tmpz * scale_p; #endif } else { /* rescale present momentum */ scale_p = sqrt( de * 2.0 * MASSE(p,i) / p_0_square + 1.0 ); IMPULS(p,i,X) *= scale_p; IMPULS(p,i,Y) *= scale_p; #ifndef TWOD IMPULS(p,i,Z) *= scale_p; #endif } } } } /* void laser_rescale_1()*/
double get_surface() { /* Calculates a 1D (anlong the x-axes) density distribution, to 'find' the surface. 1D is sufficient here, because laser radiation is only availible in (1,0,0) direction (yet) */ /* One density-cell is chosen to be ~ 2.5 A, but I depends on the system*/ double deltax = 2.5; int ndcells = (int)(box_x.x/deltax); int *xdens_1, *xdens_2 = NULL; int l,k; int rightside = ndcells; int leftside = 0; #ifdef MPI2 MPI_Alloc_mem( ndcells * sizeof(int), MPI_INFO_NULL, &xdens_1 ); MPI_Alloc_mem( ndcells * sizeof(int), MPI_INFO_NULL, &xdens_2 ); #elif defined(MPI) xdens_1 = (int *) malloc( ndcells * sizeof(int) ); xdens_2 = (int *) malloc( ndcells * sizeof(int) ); #else xdens_1 = (int *) malloc( ndcells * sizeof(int) ); xdens_2 = (int *) malloc( ndcells * sizeof(int) ); #endif for (l=0; l<(ndcells); l++) { xdens_1[l] = 0; xdens_2[l] = 0; } /* sum over the density-cells */ for (l=0; l<(ndcells); l++) { /* check if an atoms x-coordiante is in the l-th dcell */ for (k=0; k<NCELLS; k++) { cell *p; p = CELLPTR(k); int i; /* if so, add 1 atom to the xdens of the l-th dcell */ for (i=0; i < p->n; i++) { if( (ORT(p,i,X) > (double)l * deltax ) && (ORT(p,i,X) < (double)(l+1) * deltax ) ) xdens_2[l] += 1; } } } /* add the resultes for the different CPUs */ #ifdef MPI MPI_Reduce( xdens_2, xdens_1, ndcells, MPI_INT, MPI_SUM, 0, cpugrid); #else for(l=0; l<ndcells;l++) { xdens_1[l] = xdens_2[l]; } #endif if(myid==0){ /* find the right side of the sample, which is the most outer dcell with density != 0 */ for (l=ndcells-1;l>0;l--) { if( (xdens_1[l] == 0) && (xdens_1[l-1] != 0) ) { rightside = l-1; break; } } #ifdef DEBUG for (l=0; l<ndcells; l++) printf("num(%i): %i, intervall: %f - %f \n", l, xdens_1[l], l* deltax, (l+1)* deltax); #endif /* find the actual surface (not clusters or vapour) by starting from the very right side and look for two adjacent dcells with density 0 */ for (l=rightside; l>0; l--){ if((xdens_1[l]==0)&&(xdens_1[l-1]==0)) break; } leftside = l+1; /* check if there are only a few atoms inside the actual top 2 guessed 'surface dcells'; if so adjusts the surface slightly; do also for the rightside */ if ((xdens_1[leftside]<500)) { if ((xdens_1[leftside+1]<500)) leftside = l+3; else leftside = l+2; } if ((xdens_1[rightside]<500)) { if((xdens_1[rightside-1]<500)) rightside-=2; else rightside-=1; } } /* free arrays*/ #ifdef MPI2 MPI_Free_mem(xdens_1); MPI_Free_mem(xdens_2); #elif defined(MPI) free(xdens_1); free(xdens_2); #else free(xdens_1); free(xdens_2); #endif laser_atom_vol=calc_laser_atom_vol(deltax, leftside, rightside, xdens_1); /* send the new value for laser_atom_vol to the other CPUs */ #ifdef MPI MPI_Bcast( &laser_atom_vol, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD); #endif #ifdef PDECAY double samplesize = (rightside - leftside)*deltax; ramp_start = (1.0 - ramp_fraction) * samplesize + (double)(leftside+0.5)*deltax; ramp_end = (double)(rightside) * deltax; #ifdef MPI MPI_Bcast( &ramp_start, 1, REAL, 0, MPI_COMM_WORLD); MPI_Bcast( &ramp_end, 1, REAL, 0, MPI_COMM_WORLD); #endif #endif /* return the new actual surface x-coordiante (parameter laser_offset) */ return (double)(leftside+0.5)*deltax; }
void calc_forces_neb(void) { real dl2=0.0, dr2=0.0, drl=0.0, d2=0.0, f2=0.0, f2max=0.0, drlmax=0.0,df=0.0; real tmp,tmp1,tmp2, cosphi, fphi, src[3], dest[3], *d=pos; real kr,kl; int k, i; int var_k=0; int myimage,maximage; real V_previous, V_actual, V_next; real deltaVmin,deltaVmax; real normdr,normdl,inormd; real Eref,Emax,Emin,delta_E; real ratio_plus,ratio_minus,abs_next,abs_previous ; real k_sum, k_diff,tmpl,tmpr; real tmp_neb_ks[NEB_MAXNREP] INIT(zero100); real felastfact=0.0; myimage = myrank; /* get info about the energies of the different images */ neb_image_energies[ myimage]=tot_pot_energy; MPI_Allreduce(neb_image_energies , neb_epot_im, NEB_MAXNREP, REAL, MPI_SUM, MPI_COMM_WORLD); Emax=-999999999999999; Emin=999999999999999; for(i=0;i<neb_nrep;i++) { if(neb_epot_im[i]>=Emax) { Emax=neb_epot_im[i]; maximage=i; } if(neb_epot_im[i]<=Emin) { Emin=neb_epot_im[i]; } } if(steps == neb_cineb_start) { if(neb_climbing_image > 0) { if(myrank==0) { if( neb_climbing_image == maximage) printf("Starting climbing image = %d (= max_Epot = %lf)\n",neb_climbing_image, Emax); else printf("Starting climbing image = %d \n WARNING: %d != %d with max_Epot = %lf)\n", \ neb_climbing_image,neb_climbing_image,maximage, Emax); } } else { neb_climbing_image = maximage; if(myrank==0) { printf("Starting climbing image, image set to %d (= max_Epot = %lf)\n",maximage, Emax); } } } /* determine variable spring constants (jcp113 p. 9901) */ tmp_neb_ks[myimage]=0; if(myrank != 0 && myrank != neb_nrep-1) { V_previous = neb_epot_im[myimage-1]; V_actual = neb_epot_im[myimage]; V_next = neb_epot_im[myimage+1]; if ( neb_kmax > 0 & neb_kmin >0 && steps > neb_vark_start) { var_k=1; k_sum = neb_kmax + neb_kmin; k_diff = neb_kmax - neb_kmin; delta_E = Emax - Emin; if (delta_E > 1.0e-12) { tmp_neb_ks[myimage] = 0.5 *(k_sum - k_diff * cos(3.141592653589793238*( neb_epot_im[myimage] - Emin )/delta_E )); } } else { tmp_neb_ks[myimage] = neb_k; } } MPI_Allreduce(tmp_neb_ks , neb_ks, NEB_MAXNREP, REAL, MPI_SUM, MPI_COMM_WORLD); /* exchange positions with neighbor replicas */ neb_sendrecv_pos(); /* determine tangent vector and the elastic spring force */ if(myrank != 0 && myrank != neb_nrep-1) { dl2=0.0;d2=0.0;dr2=0.0; kr = 0.5 * (neb_ks[myimage]+neb_ks[myimage+1]); kl = 0.5 * (neb_ks[myimage]+neb_ks[myimage-1]); /* preparation: calculate distance to left and right immage */ for (i=0; i<DIM*natoms; i+=DIM) { vektor dr,dl; real x; dl.x = pos [i ] - pos_l[i ]; dl.y = pos [i+1] - pos_l[i+1]; dl.z = pos [i+2] - pos_l[i+2]; dr.x = pos_r[i ] - pos [i ]; dr.y = pos_r[i+1] - pos [i+1]; dr.z = pos_r[i+2] - pos [i+2]; /* apply periodic boundary conditions */ if (1==pbc_dirs.x) { x = - round( SPROD(dl,tbox_x) ); dl.x += x * box_x.x; dl.y += x * box_x.y; dl.z += x * box_x.z; x = - round( SPROD(dr,tbox_x) ); dr.x += x * box_x.x; dr.y += x * box_x.y; dr.z += x * box_x.z; } if (1==pbc_dirs.y) { x = - round( SPROD(dl,tbox_y) ); dl.x += x * box_y.x; dl.y += x * box_y.y; dl.z += x * box_y.z; x = - round( SPROD(dr,tbox_y) ); dr.x += x * box_y.x; dr.y += x * box_y.y; dr.z += x * box_y.z; } if (1==pbc_dirs.z) { x = - round( SPROD(dl,tbox_z) ); dl.x += x * box_z.x; dl.y += x * box_z.y; dl.z += x * box_z.z; x = - round( SPROD(dr,tbox_z) ); dr.x += x * box_z.x; dr.y += x * box_z.y; dr.z += x * box_z.z; } dRleft[i ] = dl.x; dRleft[i+1] = dl.y; dRleft[i+2] = dl.z; dRright[i ] = dr.x; dRright[i+1] = dr.y; dRright[i+2] = dr.z; } /* computation of the tangent requires 2 steps: determination of the direction and then normalization */ /* here we use only the improved tangent method */ if ( ( V_next > V_actual ) && ( V_actual > V_previous ) ) { for (i=0; i<DIM*natoms; i+=DIM) { tau[i ] = dRright[i ]; tau[i+1] = dRright[i+1]; tau[i+2] = dRright[i+2]; d2 += dRright[i ]*dRright[i ]; d2 += dRright[i+1]*dRright[i+1]; d2 += dRright[i+2]*dRright[i+2]; } tmp=1.0/sqrt(d2); for (i=0; i<DIM*natoms; i+=DIM) { tau[i ] *= tmp; tau[i+1] *= tmp; tau[i+2] *= tmp; if (var_k==1) { felastfact += tau[i ] * ( -kr *dRright[i ] + kl * dRleft[i ]); felastfact += tau[i+1] * ( -kr *dRright[i+1] + kl * dRleft[i+1]); felastfact += tau[i+2] * ( -kr *dRright[i+2] + kl * dRleft[i+2]); } else { felastfact += tau[i ] * ( - dRright[i ] + dRleft[i ]); felastfact += tau[i+1] * ( - dRright[i+1] + dRleft[i+1]); felastfact += tau[i+2] * ( - dRright[i+2] + dRleft[i+2]); } } } else if ( ( V_next < V_actual ) && ( V_actual < V_previous ) ) { for (i=0; i<DIM*natoms; i+=DIM) { tau[i ] = dRleft[i ]; tau[i+1] = dRleft[i+1]; tau[i+2] = dRleft[i+2]; d2 += dRleft[i ]*dRleft[i ]; d2 += dRleft[i+1]*dRleft[i+1]; d2 += dRleft[i+2]*dRleft[i+2]; } tmp=1.0/sqrt(d2); for (i=0; i<DIM*natoms; i+=DIM) { tau[i ] *= tmp; tau[i+1] *= tmp; tau[i+2] *= tmp; if (var_k==1) { felastfact += tau[i ] * ( -kr *dRright[i ] + kl * dRleft[i ]); felastfact += tau[i+1] * ( -kr *dRright[i+1] + kl * dRleft[i+1]); felastfact += tau[i+2] * ( -kr *dRright[i+2] + kl * dRleft[i+2]); } else { felastfact += tau[i ] * ( - dRright[i ] + dRleft[i ]); felastfact += tau[i+1] * ( - dRright[i+1] + dRleft[i+1]); felastfact += tau[i+2] * ( - dRright[i+2] + dRleft[i+2]); } } } else { abs_next = FABS( V_next - V_actual ); abs_previous = FABS( V_previous - V_actual ); deltaVmax = MAX( abs_next, abs_previous ); deltaVmin = MIN( abs_next, abs_previous ); for (i=0; i<DIM*natoms; i+=DIM) { dr2 += dRright[i ]*dRright[i ]; dr2 += dRright[i+1]*dRright[i+1]; dr2 += dRright[i+2]*dRright[i+2]; dl2 += dRleft[i ]*dRleft[i ]; dl2 += dRleft[i+1]*dRleft[i+1]; dl2 += dRleft[i+2]*dRleft[i+2]; } tmpl=1.0/sqrt(dl2); tmpr=1.0/sqrt(dr2); for (i=0; i<DIM*natoms; i+=DIM) { vektor dl, dr; dr.x = dRright[i ]*tmpr; dr.y = dRright[i+1]*tmpr; dr.z = dRright[i+2]*tmpr; dl.x = dRleft[i ]*tmpl; dl.y = dRleft[i+1]*tmpl; dl.z = dRleft[i+2]*tmpl; if (V_next > V_previous ) { tau[i ] = dr.x * deltaVmax + dl.x * deltaVmin ; tau[i+1] = dr.y * deltaVmax + dl.y * deltaVmin; tau[i+2] = dr.z * deltaVmax + dl.z * deltaVmin; } else if ( V_next < V_previous ) { tau[i ] = dr.x * deltaVmin + dl.x * deltaVmax ; tau[i+1] = dr.y * deltaVmin + dl.y * deltaVmax; tau[i+2] = dr.z * deltaVmin + dl.z * deltaVmax; } else { tau[i ] = dr.x + dl.x; tau[i+1] = dr.y + dl.y; tau[i+2] = dr.z + dl.z; } d2 += tau[i ]*tau[i ]; d2 += tau[i+1]*tau[i+1]; d2 += tau[i+2]*tau[i+2]; } tmp=1.0/sqrt(d2); for (i=0; i<DIM*natoms; i+=DIM) { tau[i ] *= tmp; tau[i+1] *= tmp; tau[i+2] *= tmp; if (var_k==1) { felastfact += tau[i ] * ( -kr *dRright[i ] + kl * dRleft[i ]); felastfact += tau[i+1] * ( -kr *dRright[i+1] + kl * dRleft[i+1]); felastfact += tau[i+2] * ( -kr *dRright[i+2] + kl * dRleft[i+2]); } else { felastfact += tau[i ] * ( - dRright[i ] + dRleft[i ]); felastfact += tau[i+1] * ( - dRright[i+1] + dRleft[i+1]); felastfact += tau[i+2] * ( - dRright[i+2] + dRleft[i+2]); } } } /* finally construct the spring force */ for (i=0; i<DIM*natoms; i+=DIM) { if (var_k==1) { f[i ] = - tau[i ] *felastfact; f[i+1] = - tau[i+1] *felastfact; f[i+2] = - tau[i+2] *felastfact; } else { f[i ] = - neb_k * tau[i ] *felastfact; f[i+1] = - neb_k * tau[i+1] *felastfact; f[i+2] = - neb_k * tau[i+2] *felastfact; } } }// end if(myrank != 0 && mrank != neb_nrep-1) /* calculate the neb-force */ if(myrank != 0 && myrank != neb_nrep-1) { // first scalar product of -force and tangent vector tmp = 0.0; for (k=0; k<NCELLS; k++) { cell *p = CELLPTR(k); for (i=0; i<p->n; i++) { int n = NUMMER(p,i); tmp -= tau X(n) * KRAFT(p,i,X); tmp -= tau Y(n) * KRAFT(p,i,Y); tmp -= tau Z(n) * KRAFT(p,i,Z); } } // add tmp times the tangent vector // and the spring force for (k=0; k<NCELLS; k++) { cell *p = CELLPTR(k); for (i=0; i<p->n; i++) { int n = NUMMER(p,i); if(myimage == neb_climbing_image && (steps >= neb_cineb_start)) { KRAFT(p,i,X) += 2.0*tmp * tau X(n); KRAFT(p,i,Y) += 2.0*tmp * tau Y(n); KRAFT(p,i,Z) += 2.0*tmp * tau Z(n); } else { KRAFT(p,i,X) += tmp * tau X(n) + f X(n); KRAFT(p,i,Y) += tmp * tau Y(n) + f Y(n); KRAFT(p,i,Z) += tmp * tau Z(n) + f Z(n); } } } } // if(myrank != 0 && mrank != neb_nrep-1) }
void init_fcs(void) { FCSResult res; fcs_int srf = 1; char *method; fcs_int pbc [3] = { pbc_dirs.x, pbc_dirs.y, pbc_dirs.z }; fcs_float BoxX[3] = { box_x.x, box_x.y, box_x.z }; fcs_float BoxY[3] = { box_y.x, box_y.y, box_y.z }; fcs_float BoxZ[3] = { box_z.x, box_z.y, box_z.z }; fcs_float off [3] = { 0.0, 0.0, 0.0 }; /* subtract CM momentum */ if (0 == imdrestart) { int i, k; real ptot[4], ptot_2[4], px, py, pz; ptot[0] = 0.0; ptot[1] = 0.0; ptot[2] = 0.0, ptot[3] = 0.0; for (k=0; k<NCELLS; ++k) { /* loop over all cells */ cell *p = CELLPTR(k); for (i=0; i<p->n; i++) { ptot[0] += IMPULS(p,i,X); ptot[1] += IMPULS(p,i,Y); ptot[2] += IMPULS(p,i,Z); ptot[3] += MASSE(p,i); } } #ifdef MPI MPI_Allreduce( ptot, ptot_2, 4, REAL, MPI_SUM, cpugrid); ptot[0] = ptot_2[0]; ptot[1] = ptot_2[1]; ptot[2] = ptot_2[2]; ptot[3] = ptot_2[3]; #endif px = ptot[0]/ptot[3]; py = ptot[1]/ptot[3]; pz = ptot[2]/ptot[3]; for (k=0; k<NCELLS; ++k) { /* loop over all cells */ cell *p = CELLPTR(k); for (i=0; i<p->n; i++) { IMPULS(p,i,X) -= px * MASSE(p,i); IMPULS(p,i,Y) -= py * MASSE(p,i); IMPULS(p,i,Z) -= pz * MASSE(p,i); } } } switch (fcs_method) { case FCS_METH_DIRECT: method = "direct"; break; case FCS_METH_PEPC: method = "pepc"; break; case FCS_METH_FMM: method = "fmm"; break; case FCS_METH_P3M: method = "p3m"; srf = fcs_near_field_flag; break; case FCS_METH_P2NFFT: method = "p2nfft"; srf = fcs_near_field_flag; break; case FCS_METH_VMG: method = "vmg"; break; case FCS_METH_PP3MG: method = "pp3mg"; break; } /* initialize handle and set common parameters */ res = fcs_init(&handle, method, cpugrid); ASSERT_FCS(res); res = fcs_set_common(handle, srf, BoxX, BoxY, BoxZ, off, pbc, natoms); ASSERT_FCS(res); res = fcs_require_virial(handle, 1); ASSERT_FCS(res); /* set method specific parameters */ switch (fcs_method) { #ifdef FCS_ENABLE_DIRECT case FCS_METH_DIRECT: /* nothing to do */ break; #endif #ifdef FCS_ENABLE_PEPC case FCS_METH_PEPC: res = fcs_pepc_setup(handle, (fcs_float)fcs_pepc_eps, (fcs_float)fcs_pepc_theta ); ASSERT_FCS(res); res = fcs_pepc_set_num_walk_threads( handle, (fcs_int)fcs_pepc_nthreads ); ASSERT_FCS(res); break; #endif #ifdef FCS_ENABLE_FMM case FCS_METH_FMM: res = fcs_fmm_set_absrel(handle, (fcs_int)fcs_fmm_absrel); ASSERT_FCS(res); res = fcs_fmm_set_tolerance_energy(handle, (fcs_float)fcs_tolerance); ASSERT_FCS(res); break; #endif #ifdef FCS_ENABLE_P3M case FCS_METH_P3M: if (0==srf) { res = fcs_p3m_set_r_cut(handle, (fcs_float)fcs_rcut); ASSERT_FCS(res); } res = fcs_set_tolerance(handle, FCS_TOLERANCE_TYPE_FIELD, (fcs_float)fcs_tolerance); ASSERT_FCS(res); if (fcs_grid_dim.x) res = fcs_p3m_set_grid(handle, (fcs_int)fcs_grid_dim.x); ASSERT_FCS(res); break; #endif #ifdef FCS_ENABLE_P2NFFT case FCS_METH_P2NFFT: if (0==srf) { res = fcs_p2nfft_set_r_cut(handle, (fcs_float)fcs_rcut); ASSERT_FCS(res); } res = fcs_set_tolerance(handle, FCS_TOLERANCE_TYPE_FIELD, (fcs_float)fcs_tolerance); ASSERT_FCS(res); if (fcs_grid_dim.x) { res = fcs_p2nfft_set_grid(handle, (fcs_int)fcs_grid_dim.x, (fcs_int)fcs_grid_dim.y, (fcs_int)fcs_grid_dim.z); ASSERT_FCS(res); } if (fcs_p2nfft_intpol_order) { res = fcs_p2nfft_set_pnfft_interpolation_order(handle, (fcs_int)fcs_p2nfft_intpol_order); ASSERT_FCS(res); } if (fcs_p2nfft_epsI) { res = fcs_p2nfft_set_epsI(handle, (fcs_float)fcs_p2nfft_epsI); ASSERT_FCS(res); } //res = fcs_p2nfft_set_pnfft_window_by_name(handle, "bspline"); //ASSERT_FCS(res); break; #endif #ifdef FCS_ENABLE_VMG case FCS_METH_VMG: if (fcs_vmg_near_field_cells) { res = fcs_vmg_set_near_field_cells(handle, (fcs_int)fcs_vmg_near_field_cells); ASSERT_FCS(res); } if (fcs_vmg_interpol_order) { res = fcs_vmg_set_interpolation_order(handle, (fcs_int)fcs_vmg_interpol_order); ASSERT_FCS(res); } if (fcs_vmg_discr_order) { res = fcs_vmg_set_discretization_order(handle, (fcs_int)fcs_vmg_discr_order); ASSERT_FCS(res); } if (fcs_iter_tolerance > 0) { res = fcs_vmg_set_precision(handle, (fcs_float)fcs_iter_tolerance); ASSERT_FCS(res); } break; #endif #ifdef FCS_ENABLE_PP3MG case FCS_METH_PP3MG: /* use default values, if not specified otherwise */ if (fcs_grid_dim.x) { res = fcs_pp3mg_set_cells_x(handle, (fcs_int)fcs_grid_dim.x); ASSERT_FCS(res); res = fcs_pp3mg_set_cells_y(handle, (fcs_int)fcs_grid_dim.y); ASSERT_FCS(res); res = fcs_pp3mg_set_cells_z(handle, (fcs_int)fcs_grid_dim.z); ASSERT_FCS(res); } if (fcs_pp3mg_ghosts) { res = fcs_pp3mg_set_ghosts(handle, (fcs_int)fcs_pp3mg_ghosts); ASSERT_FCS(res); } if (fcs_pp3mg_degree) { res = fcs_pp3mg_set_degree(handle, (fcs_int)fcs_pp3mg_degree); ASSERT_FCS(res); } if (fcs_pp3mg_max_part) { res = fcs_pp3mg_set_max_particles(handle, (fcs_int)fcs_pp3mg_max_part); ASSERT_FCS(res); } if (fcs_max_iter) { res = fcs_pp3mg_set_max_iterations(handle,(fcs_int)fcs_max_iter); ASSERT_FCS(res); } if (fcs_iter_tolerance > 0) { res = fcs_pp3mg_set_tol(handle, (fcs_float)fcs_iter_tolerance); ASSERT_FCS(res); } break; #endif default: error("FCS method unknown or not implemented"); break; } pack_fcs(); res = fcs_tune(handle, nloc, nloc_max, pos, chg); ASSERT_FCS(res); /* inform about tuned parameters */ switch (fcs_method) { fcs_int grid_dim[3]; fcs_float r_cut; #ifdef FCS_ENABLE_P3M case FCS_METH_P3M: res = fcs_p3m_get_r_cut(handle, &r_cut); ASSERT_FCS(res); res = fcs_p3m_get_grid(handle, grid_dim); ASSERT_FCS(res); if (0==myid) printf("FCS: Tuned grid dimensions, cutoff: %d %d %d, %f\n", grid_dim[0], grid_dim[1], grid_dim[2], r_cut); break; #endif #ifdef FCS_ENABLE_P2NFFT case FCS_METH_P2NFFT: res = fcs_p2nfft_get_grid(handle, grid_dim, grid_dim+1, grid_dim+2); ASSERT_FCS(res); res = fcs_p2nfft_get_r_cut(handle, &r_cut); ASSERT_FCS(res); if (0==myid) printf("FCS: Tuned grid dimensions, cutoff: %d %d %d, %f\n", grid_dim[0], grid_dim[1], grid_dim[2], r_cut); break; #endif #ifdef FCS_ENABLE_PP3MG case FCS_METH_PP3MG: res = fcs_pp3mg_get_cells_x(handle, grid_dim ); ASSERT_FCS(res); res = fcs_pp3mg_get_cells_y(handle, grid_dim+1); ASSERT_FCS(res); res = fcs_pp3mg_get_cells_z(handle, grid_dim+2); if (0==myid) printf("FCS: Tuned grid dimensions: %d %d %d\n", grid_dim[0], grid_dim[1], grid_dim[2]); ASSERT_FCS(res); break; #endif default: break; } /* add near-field potential, after fcs_tune */ if (0==srf) fcs_update_pottab(); }
void unpack_fcs(void) { fcs_float vir[9] = { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 }; FCSResult result; real pot1, pot2, e, c, sum=0.0, fac=0.5; int n, m, k, i; /* extract output and distribute it to cell array */ n=0; m=0; pot1=0.0; for (k=0; k<NCELLS; ++k) { cell *p = CELLPTR(k); for (i=0; i<p->n; ++i) { c = CHARGE(p,i) * coul_eng; KRAFT(p,i,X) += field[n++] * c; KRAFT(p,i,Y) += field[n++] * c; KRAFT(p,i,Z) += field[n++] * c; e = pot[m++] * c * fac; POTENG(p,i) += e; pot1 += e; } } /* unpack virial */ result = fcs_get_virial(handle, vir); ASSERT_FCS(result); #ifdef P_AXIAL vir_xx += vir[0]; vir_yy += vir[4]; vir_zz += vir[8]; #else virial += vir[0] + vir[4] + vir[8]; #endif #ifdef STRESS_TENS if (do_press_calc) { /* distribute virial tensor evenly on all atoms */ sym_tensor pp; pp.xx = vir[0] / natoms; pp.yy = vir[4] / natoms; pp.zz = vir[8] / natoms; pp.yz = (vir[5]+vir[7]) / (2*natoms); pp.zx = (vir[2]+vir[6]) / (2*natoms); pp.xy = (vir[1]+vir[3]) / (2*natoms); for (k=0; k<NCELLS; ++k) { cell *p = CELLPTR(k); for (i=0; i<p->n; ++i) { PRESSTENS(p,i,xx) += pp.xx; PRESSTENS(p,i,yy) += pp.yy; PRESSTENS(p,i,zz) += pp.zz; PRESSTENS(p,i,yz) += pp.yz; PRESSTENS(p,i,zx) += pp.zx; PRESSTENS(p,i,xy) += pp.xy; } } } #endif /* sum up potential energy */ #ifdef MPI MPI_Allreduce( &pot1, &pot2, 1, MPI_DOUBLE, MPI_SUM, cpugrid); tot_pot_energy += pot2; #else tot_pot_energy += pot1; #endif }