/* Return the molecule orbital values */ void Average_ekt::calc_mos(Sample_point * sample, int e, Array2 <dcomplex> & movals) { if(complex_orbitals) { movals.Resize(nmo,1); cmomat->updateVal(sample,e,0,movals); } else { Array2 <doublevar> movals_d(nmo,1); momat->updateVal(sample,e,0,movals_d); movals.Resize(nmo,1); for(int i=0; i< nmo; i++) { movals(i,0)=movals_d(i,0); } } }
/* Calculate the values, gradient and laplacians of molecule orbitals, returned as: [val, grad, lap] */ void Average_ekt::calc_mosLap(Sample_point * sample, int e, Array2 <dcomplex> & molaps) { if(complex_orbitals) { molaps.Resize(nmo,5); cmomat->updateLap(sample, e, 0, molaps); } else { Array2 <doublevar> molaps_d(nmo,5); momat->updateLap(sample,e,0,molaps_d); molaps.Resize(nmo,5); for(int i=0; i< nmo; i++) { for (int d=0; d<5; d++) molaps(i,d)=molaps_d(i,d); } } }
void Minimization_function::init_overlap_ma_R(Array2<doublevar> &overlap_ma_R){ assert(overlap_ma_R.GetDim(0)==overlap_ma_R.GetDim(1)); _overlap_ma_R.Resize(overlap_ma_R.GetDim(0),overlap_ma_R.GetDim(0)); for(int k=0;k<overlap_ma_R.GetDim(0);k++) for(int l=0;l<=k;l++) _overlap_ma_R(l,k)=_overlap_ma_R(k,l)=overlap_ma_R(k,l); }
int Molecular_system::getBounds(Array2 <doublevar> & latticevec,Array1<doublevar> & origin) { cout << "getBounds()" << endl; latticevec.Resize(3,3); latticevec=0.0; Array1 <doublevar> minmax(6); for(int d=0; d< 3; d++) { minmax(2*d)=minmax(2*d+1)=0.0; } int nions=nIons(); Array1 <doublevar> ionpos(3); for(int i=0; i< nions; i++) { getIonPos(i,ionpos); for(int d=0; d< 3; d++) { if(ionpos(d) < minmax(2*d)) minmax(2*d) = ionpos(d); if(ionpos(d) > minmax(2*d+1)) minmax(2*d+1)=ionpos(d); } } for(int d=0; d< 3; d++) { // minmax(2*d)-=4.0; //minmax(2*d+1)+=4.0; origin(d)=minmax(2*d)-4.0; latticevec(d,d)=minmax(2*d+1)-origin(d)+4.0; } cout << "getBounds done" << endl; return 1; }
//---------------------------------------------------------------------- //Keeping this separate from calcLoc for efficiency reasons.. //It's most likely faster to add to a double than fill matrices.. //We also don't need to calculate the ion-ion interaction for this void Molecular_system::separatedLocal(Sample_point * sample,Array1 <doublevar> & onebody, Array2 <doublevar> & twobody) { int nions=sample->ionSize(); int nelectrons=sample->electronSize(); onebody.Resize(nelectrons); twobody.Resize(nelectrons,nelectrons); onebody=0.0; twobody=0.0; Array1 <doublevar> R(5); sample->updateEIDist(); sample->updateEEDist(); for(int e=0; e< nelectrons; e++) { for(int i=0; i < nions; i++) { sample->getEIDist(e,i, R); onebody(e)-=sample->getIonCharge(i)/R(0); } } Array1 <doublevar> R2(5); for(int i=0; i< nelectrons; i++) { for(int j=i+1; j<nelectrons; j++) { sample->getEEDist(i,j,R2); twobody(i,j)+= 1/R2(0); } } Array1 <doublevar> pos(3); for(int e=0; e< nelectrons; e++) { sample->getElectronPos(e,pos); for(int d=0; d< 3; d++) onebody(e)-=electric_field(d)*pos(d); } }
void GeneralizedEigenSystemSolverRealGeneralMatrices(Array2 < doublevar > & Ain, Array1 <dcomplex> & W, Array2 <doublevar> & VL, Array2 <doublevar> & VR) { #ifdef USE_LAPACK //if LAPACK int N=Ain.dim[0]; Array2 <doublevar> A_temp=Ain; //,VL(N,N),VR(N,N); Array1 <doublevar> WORK,RWORK(2*N),WI(N),WR(N); WI.Resize(N); VL.Resize(N,N); VR.Resize(N,N); int info; int NB=64; int NMAX=N; int lda=NMAX; int ldb=NMAX; int LWORK=5*NMAX; WORK.Resize(LWORK); info= dgeev('V','V',N,A_temp.v, lda,WR.v,WI.v,VL.v,lda,VR.v,lda,WORK.v,LWORK); if(info>0) error("Internal error in the LAPACK routine dgeev",info); if(info<0) error("Problem with the input parameter of LAPACK routine dgeev in position ",-info); W.Resize(N); for(int i=0; i< N; i++) { W(i)=dcomplex(WR(i),WI(i)); } // for (int i=0; i<N; i++) // evals(i)=W[N-1-i]; // for (int i=0; i<N; i++) { // for (int j=0; j<N; j++) { // evecs(j,i)=A_temp(N-1-i,j); // } // } //END OF LAPACK #else //IF NO LAPACK error("need LAPACK for eigensystem solver for general matrices"); #endif //END OF NO LAPACK }
void Properties_final_average::twoPointForces(Array2 <doublevar> & forces ) { int nwf=avg.GetDim(1); int naux=aux_energy.GetDim(0); if(naux%2!=0) error("there can't be two-point forces, since" " there's an odd number of differences"); forces.Resize(naux/2, 2); for(int i=0; i < naux; i+=2) { assert(aux_size(i)==aux_size(i+1)); forces(i/2,0)=(aux_diff(i,0)-aux_diff(i+1,0))/(aux_size(i)*2); forces(i/2,1)=sqrt( (aux_differr(i,0)+aux_differr(i+1,0))/4)/aux_size(i); } }
//---------------------------------------------------------------------- // // //Note this needs to be changed for non-zero k-points! doublevar Average_ekt::gen_sample(int nstep, doublevar tstep, int e, Array2 <dcomplex> & movals, Sample_point * sample) { int ndim=3; Array1 <doublevar> r(ndim),rold(ndim); Array2 <dcomplex> movals_old(nmo,1); movals.Resize(nmo,1); sample->getElectronPos(e,rold); calc_mos(sample,e,movals_old); doublevar acc=0; for(int step=0; step < nstep; step++) { for(int d=0; d< ndim; d++) { r(d)=rold(d)+sqrt(tstep)*rng.gasdev(); } sample->setElectronPos(e,r); calc_mos(sample,e,movals); doublevar sum_old=0,sum=0; for(int mo=0; mo < nmo; mo++) { sum_old+=norm(movals_old(mo,0)); sum+=norm(movals(mo,0)); } if(rng.ulec() < sum/sum_old) { movals_old=movals; rold=r; acc++; } } //cout << "acceptance " << acc/nstep << endl; movals=movals_old; sample->setElectronPos(e,rold); doublevar sum=0; for(int mo=0; mo < nmo; mo++) { sum+=norm(movals(mo,0)); } return sum; }
/*! Uses flat file of form: ncenters Label x y z Label x y z ... */ void read_centerpos(string & filename, Array2 <doublevar> & position, vector <string> & labels) { int ncenters; ifstream centin(filename.c_str()); if(!centin) error("Couldn't open ", filename); centin >> ncenters; labels.clear(); position.Resize(ncenters,3); string labeltemp; for(int i=0; i< ncenters; i++) { centin >> labeltemp; labels.push_back(labeltemp); for(int d=0; d< 3; d++) { centin >> position(i,d); } } centin.close(); }
void get_onebody_parms_diffspin(vector<string> & words, vector<string> & atomnames, Array2 <doublevar> & unique_parameters, int s, Array1 <int> & nparms, vector <string> & parm_labels, Array2 <int> & linear_parms, int & counter, Array1 <int> & parm_centers) { vector < vector < string> > parmtxt; vector <string> parmtmp; unsigned int pos=0; if(s==0){ while(readsection(words, pos,parmtmp, "LIKE_COEFFICIENTS")) parmtxt.push_back(parmtmp); if(parmtxt.size()==0) error("Didn't find LIKE COEFFICIENTS in one of threebody spin Jastrow section"); } else{ while(readsection(words, pos,parmtmp, "UNLIKE_COEFFICIENTS")) parmtxt.push_back(parmtmp); if(parmtxt.size()==0) error("Didn't find UNLIKE COEFFICIENTS in one of threebody spin Jastrow section"); } int nsec=parmtxt.size(); int maxsize=0; nparms.Resize(nsec); for(int i=0; i< nsec; i++) nparms(i)=parmtxt[i].size()-1; for(int i=0; i< nsec; i++) if(maxsize < nparms(i)) maxsize=nparms(i); unique_parameters.Resize(nsec, maxsize); linear_parms.Resize(nsec, maxsize); //cout <<"spin "<<s<<" counter "<<counter<<endl; for(int i=0; i< nsec; i++) { for(int j=0; j< nparms(i); j++) { unique_parameters(i,j)=atof(parmtxt[i][j+1].c_str()); linear_parms(i,j)=counter++; } } parm_labels.resize(nsec); for(int i=0; i< nsec; i++) parm_labels[i]=parmtxt[i][0]; int natoms=atomnames.size(); parm_centers.Resize(natoms); parm_centers=-1; for(int i=0; i< nsec; i++) { int found=0; for(int j=0; j< natoms; j++) { if(atomnames[j]==parm_labels[i]) { parm_centers(j)=i; found=1; } } if(!found) error("Couldn't find a matching center for ", parm_labels[i]); } doublevar scale; if(readvalue(words, pos=0, scale, "SCALE")) { for(int i=0; i< unique_parameters.GetDim(0); i++) { for(int j=0; j< nparms(i); j++) { unique_parameters(i,j)*=scale; } } } }
//A is upper triangular and gives the angle to rotate for a pair i,j void make_rotation_matrix(const Array2 <doublevar> & A, Array2 <doublevar> & B) { int n=A.GetDim(0); assert(A.GetDim(1)==n); B.Resize(n,n); B=0.; for(int i=0; i< n; i++) B(i,i)=1.0; Array2 <doublevar> tmp(n,n),tmp2(n,n); Array1 <doublevar> tmpveci(n),tmpvecj(n); doublevar s,co; for(int i=0; i< n; i++) { for(int j=i+1; j < n; j++) { if(fabs(A(i,j)) > 1e-12) { co=cos(A(i,j)); s=sin(A(i,j)); /* tmp2=0.0; tmp=0.0; for(int k=0; k< n; k++) tmp2(k,k)=1.0; tmp2(i,i)=co; tmp2(i,j)=-s; tmp2(j,i)=s; tmp2(j,j)=co; //very very inefficient, for testing.. //MultiplyMatrices(tmp2,B,tmp,n); tmp=0; for(int a=0; a< n; a++) { for(int c=0; c< n; c++) { tmp(a,c)+=tmp2(a,a)*B(a,c); } } for(int c=0; c< n; c++) { tmp(i,c)+=tmp2(i,j)*B(j,c); tmp(j,c)+=tmp2(j,i)*B(i,c); } B=tmp; */ //------ for(int c=0; c< n; c++) { tmpveci(c)=B(i,c); tmpvecj(c)=B(j,c); } for(int c=0; c < n; c++) { B(i,c)*=co; B(j,c)*=co; } for(int c=0; c< n; c++) { B(i,c)+= -s*tmpvecj(c); B(j,c)+= s*tmpveci(c); } } } } /* cout << "rotation matrix " << endl; for(int i=0; i< n; i++) { for(int j=0; j< n; j++) { cout << B(i,j) << " "; } cout << endl; } */ }
//---------------------------------------------------------------------- void Wannier_method::optimize_rotation(Array3 <dcomplex> & eikr, Array2 <doublevar> & R ) { int norb=eikr.GetDim(1); Array2 <doublevar> gvec(3,3); sys->getPrimRecipLattice(gvec); Array1 <doublevar> gnorm(3); gnorm=0; for(int d=0; d < 3; d++) { for(int d1=0; d1 < 3; d1++) { gnorm(d)+=gvec(d,d1)*gvec(d,d1); } gnorm(d)=sqrt(gnorm(d)); } for(int i=0; i< norb; i++) { cout << "rloc2 " << i << " "; for(int d=0; d< 3; d++) { cout << -log(norm(eikr(d,i,i)))/(gnorm(d)*gnorm(d)) << " "; } cout << endl; } Array2 <doublevar> Rgen(norb,norb),Rgen_save(norb,norb); //R(norb,norb); R.Resize(norb,norb); //Array2 <dcomplex> tmp(norb,norb),tmp2(norb,norb); //Shake up the angles, since often the original orbitals //are at a maximum and derivatives are zero. Array2 <doublevar> deriv(norb,norb); Rgen=0.0; for(int ii=0; ii< norb; ii++) { for(int jj=ii+1; jj< norb; jj++) { Rgen(ii,jj)=rng.gasdev()*pi*shake; } } for(int step=0; step < max_opt_steps; step++) { doublevar fbase=evaluate_local(eikr,Rgen,R); for(int ii=0; ii <norb; ii++) { cout << "deriv "; for(int jj=ii+1; jj < norb; jj++) { doublevar save_rgeniijj=Rgen(ii,jj); doublevar h=1e-6; Rgen(ii,jj)+=h; doublevar func=evaluate_local(eikr,Rgen,R); deriv(ii,jj)=(func-fbase)/h; Rgen(ii,jj)=save_rgeniijj; cout << deriv(ii,jj) << " "; } cout << endl; } doublevar rloc_thresh=0.0001; Rgen_save=Rgen; doublevar best_func=1e99, best_tstep=0.0; doublevar bracket_tstep=0.0; doublevar last_func=fbase; for(doublevar tstep=0.01; tstep < 20.0; tstep*=2.0) { doublevar func=eval_tstep(eikr,Rgen,Rgen_save,deriv,tstep,R); cout << "tstep " << tstep << " func " << func << endl; if(func > fbase or func > last_func) { bracket_tstep=tstep; break; } else last_func=func; } cout << "bracket_tstep " << bracket_tstep << endl; doublevar resphi=2.-(1.+sqrt(5.))/2.; doublevar a=0, b=resphi*bracket_tstep, c=bracket_tstep; doublevar af=fbase, bf=eval_tstep(eikr,Rgen,Rgen_save,deriv,b,R), cf=eval_tstep(eikr,Rgen,Rgen_save,deriv,bracket_tstep,R); cout << "first step a,b,c " << a << " " << b << " " << c << " funcs " << af << " " << bf << " " << cf << endl; for(int it=0; it < 20; it++) { doublevar d,df; if( (c-b) > (b-a)) d=b+resphi*(c-b); else d=b-resphi*(b-a); df=eval_tstep(eikr,Rgen,Rgen_save,deriv,d,R); if(df < bf) { if( (c-b) > (b-a) ) { a=b; af=bf; b=d; bf=df; } else { c=b; cf=bf; b=d; bf=df; } } else { if( (c-b) > (b-a) ) { c=d; cf=df; } else { a=d; af=df; } } cout << "step " << it << " a,b,c " << a << " " << b << " " << c << " funcs " << af << " " << bf << " " << cf << endl; } best_tstep=b; /* bool made_move=false; while (!made_move) { for(doublevar tstep=0.00; tstep < max_tstep; tstep+=0.1*max_tstep) { for(int ii=0; ii< norb;ii++) { for(int jj=ii+1; jj < norb; jj++) { Rgen(ii,jj)=Rgen_save(ii,jj)-tstep*deriv(ii,jj); } } doublevar func=evaluate_local(eikr,Rgen,R); if(func < best_func) { best_func=func; best_tstep=tstep; } cout << " tstep " << tstep << " " << func << endl; } if(abs(best_tstep) < 0.2*max_tstep) max_tstep*=0.5; else if(abs(best_tstep-max_tstep) < 1e-14) max_tstep*=2.0; else made_move=true; } */ for(int ii=0; ii< norb; ii++) { for(int jj=ii+1; jj < norb; jj++) { Rgen(ii,jj)=Rgen_save(ii,jj)-best_tstep*deriv(ii,jj); } } doublevar func2=evaluate_local(eikr,Rgen,R); doublevar max_change=0; for(int ii=0; ii < norb; ii++) { for(int jj=ii+1; jj< norb; jj++) { doublevar change=abs(Rgen(ii,jj)-Rgen_save(ii,jj)); if(change > max_change) max_change=change; } } cout << "tstep " << best_tstep << " rms " << sqrt(func2) << " bohr max change " << max_change <<endl; doublevar threshold=0.0001; if(max_change < threshold) break; if(abs(best_func-fbase) < rloc_thresh) break; /* bool moved=false; for(int ii=0; ii< norb; ii++) { for(int jj=ii+1; jj< norb; jj++) { doublevar save_rgeniijj=Rgen(ii,jj); doublevar best_del=0; doublevar best_f=1e99; for(doublevar del=-0.5; del < 0.5; del+=0.05) { cout << "############ for del = " << del << endl; Rgen(ii,jj)=save_rgeniijj+del; doublevar func=evaluate_local(eikr,Rgen,R); if(func < best_f) { best_f=func; best_del=del; } cout << "func " << func << endl; } Rgen(ii,jj)=save_rgeniijj+best_del; if(abs(best_del) > 1e-12) moved=true; } } if(!moved) break; */ } make_rotation_matrix(Rgen,R); }
//---------------------------------------------------------------------- void Wannier_method::calculate_overlap(Array1 <int> & orb_list, Array3 <dcomplex> & eikr, Array2 <doublevar> & phi2phi2) { Array1 <Array1 <int> > tmp_list(1); tmp_list(0)=orb_list; mymomat->buildLists(tmp_list); int list=0; Array2 <doublevar> gvec(3,3); if(!sys->getPrimRecipLattice(gvec) ) error("Need getPrimRecipLattice"); int norb=orb_list.GetDim(0); eikr.Resize(3,norb,norb); phi2phi2.Resize(norb,norb); eikr=0.0; phi2phi2=0.0; Array1 <doublevar> prop(3); Array1 <doublevar> n_lat(3); n_lat=0.0; for(int d=0; d < 3; d++) { for(int d1=0; d1 < 3; d1++) n_lat(d)+=LatticeVec(d,d1)*LatticeVec(d,d1); n_lat(d)=sqrt(n_lat(d)); prop(d)=resolution/double(n_lat(d)); //cout << "prop " << prop(d) << " res " << resolution << " norm " << n_lat(d) << endl; } cout << "calculating...." << endl; mymovals.Resize(norb,1); Array1 <doublevar> r(3),x(3); Array1 <doublevar> norm_orb(norb); norm_orb=0.0; int totpts=0; for(r(0)=origin(0); r(0) < 1.0; r(0)+=prop(0)) { cout << "percent done: " << r(0)*100 << endl; for(r(1)=origin(1); r(1) < 1.0; r(1)+=prop(1)) { for(r(2)=origin(2); r(2) < 1.0; r(2)+=prop(2)) { totpts++; for(int d=0; d< 3; d++) { x(d)=0.0; for(int d1=0; d1 < 3; d1++) { x(d)+=r(d1)*LatticeVec(d1,d); } } mywalker->setElectronPos(0,x); mymomat->updateVal(mywalker,0,list,mymovals); for(int i=0; i < norb; i++) norm_orb(i)+=mymovals(i,0)*mymovals(i,0); for(int d=0; d< 3; d++) { doublevar gdotr=0; for(int d1=0; d1 < 3; d1++) { gdotr+=gvec(d,d1)*x(d1); } dcomplex exp_tmp=exp(dcomplex(0,-1)*2.*pi*gdotr); for(int i=0; i< norb; i++) { for(int j=0; j< norb; j++) { eikr(d,i,j)+=exp_tmp*mymovals(i,0)*mymovals(j,0); } } } for(int i=0; i< norb; i++) { for(int j=0; j< norb; j++) { phi2phi2(i,j)+=mymovals(i,0)*mymovals(i,0)*mymovals(j,0)*mymovals(j,0); } } } } } for(int i=0; i< norb; i++) norm_orb(i)=sqrt(norm_orb(i)); for(int d=0; d< 3; d++) { for(int i=0; i < norb; i++) { for(int j=0; j< norb; j++) { eikr(d,i,j)/=norm_orb(i)*norm_orb(j); //cout << d << " " << i << " " << j << " " << eikr(d,i,j) << endl; } } } cout << "square overlap " << endl; for(int i=0; i < norb; i++) { for(int j=0; j< norb; j++) { phi2phi2(i,j)/=(norm_orb(i)*norm_orb(i)*norm_orb(j)*norm_orb(j)); cout << sqrt(phi2phi2(i,j)) << " "; } cout << endl; } }
//A is an antisymmetric matrix and B is the output rotation matrix void make_rotation_matrix_notworking(const Array2 <doublevar> & A, Array2 <doublevar> & B) { int n=A.GetDim(0); assert(A.GetDim(1)==n); B.Resize(n,n); Array2 <dcomplex> skew(n,n),VL(n,n),VR(n,n); Array1 <dcomplex> evals(n); for(int i=0; i< n; i++) { for(int j=0; j < n; j++) { skew(i,j)=A(i,j); } } GeneralizedEigenSystemSolverComplexGeneralMatrices(skew,evals,VL,VR); cout << "evals " << endl; for(int i=0; i< n; i++) cout << evals(i) << " "; cout << endl; cout << "VR " << endl; for(int i=0; i< n; i++) { for(int j=0; j< n; j++) { cout << VR(i,j) << " "; } cout << endl; } cout << "VL " << endl; for(int i=0; i< n; i++) { for(int j=0; j< n; j++) { cout << VL(i,j) << " "; } cout << endl; } //this is horribly inefficient,most likely skew=dcomplex(0.0,0.); //we don't need that any more so we reuse it Array2 <dcomplex> work(n,n); work=dcomplex(0.0,0.); for(int i=0; i< n; i++) { skew(i,i)=exp(evals(i)); } for(int i=0; i< n; i++) { for(int j=0; j<n; j++) { for(int k=0; k< n; k++) { work(i,k)+=skew(i,j)*VR(j,k); } } } skew=dcomplex(0.,0.); for(int i=0; i< n; i++) { for(int j=0; j<n; j++) { for(int k=0; k< n; k++) { //skew(i,k)+=conj(VL(i,j))*work(j,k); skew(i,k)=conj(VR(j,i))*work(j,k); } } } cout << "rotation " << endl; for(int i=0; i< n; i++) { for(int j=0; j< n; j++) { cout << skew(i,j) << " "; } cout << endl; } }
void MO_matrix_blas::updateLap( Sample_point * sample, int e, int listnum, Array2 <doublevar> & newvals) { #ifdef USE_BLAS int centermax=centers.size(); assert(e < sample->electronSize()); assert(newvals.GetDim(1) >=5); Array1 <doublevar> R(5); static Array2 <doublevar> symmvals_temp(maxbasis,5); static Array2 <doublevar> newvals_T; Array2 <doublevar> & moCoefftmp(moCoeff_list(listnum)); int totbasis=moCoefftmp.GetDim(0); int nmo_list=moCoefftmp.GetDim(1); newvals_T.Resize(5, nmo_list); newvals_T=0.0; int b; int totfunc=0; centers.updateDistance(e, sample); static Array1 <MOBLAS_CalcObjLap> calcobjs(totbasis); int ncalcobj=0; for(int ion=0; ion < centermax; ion++) { centers.getDistance(e, ion, R); for(int n=0; n< centers.nbasis(ion); n++) { b=centers.basis(ion, n); if(R(0) < obj_cutoff(b)) { basis(b)->calcLap(R, symmvals_temp); int imax=nfunctions(b); for(int i=0; i< imax; i++) { if(R(0) < cutoff(totfunc)) { calcobjs(ncalcobj).moplace=moCoefftmp.v+totfunc*nmo_list; for(int j=0; j< 5; j++) { calcobjs(ncalcobj).sval[j]=symmvals_temp(i,j); //cblas_daxpy(nmo_list,symmvals_temp(i,j), // moCoefftmp.v+totfunc*nmo_list,1, // newvals_T.v+j*nmo_list,1); } ncalcobj++; } totfunc++; } } else { totfunc+=nfunctions(b); } } } for(int i=0; i< ncalcobj; i++) { for(int j=0; j< 5; j++) { cblas_daxpy(nmo_list,calcobjs(i).sval[j],calcobjs(i).moplace,1, newvals_T.v+j*nmo_list,1); } //cblas_dger(CblasRowMajor,5,nmo_list,1., // calcobjs(i).sval,1, // calcobjs(i).moplace,1, // newvals_T.v,nmo_list); } for(int m=0; m < nmo_list; m++) { for(int j=0; j< 5; j++) { newvals(m,j)=newvals_T(j,m); } } #else error("BLAS libraries not compiled in, so you cannot use BLAS_MO"); #endif }