@@ -94,6 +94,7 @@ void ESolver_SDFT_PW::check_che(const int nche_in)
9494void ESolver_SDFT_PW::sKG (const int nche_KG, const double fwhmin, const double wcut,
9595 const double dw_in, const int times)
9696{
97+ ModuleBase::TITLE (this ->classname ," sKG" );
9798 ModuleBase::timer::tick (this ->classname ," sKG" );
9899 cout<<" Calculating conductivity...." <<endl;
99100
@@ -456,6 +457,8 @@ void ESolver_SDFT_PW::sKG(const int nche_KG, const double fwhmin, const double w
456457
457458void ESolver_SDFT_PW:: caldos( const int nche_dos, const double sigmain, const double emin, const double emax, const double de, const int npart)
458459{
460+ ModuleBase::TITLE (this ->classname ," caldos" );
461+ ModuleBase::timer::tick (this ->classname ," caldos" );
459462 cout<<" =========================" <<endl;
460463 cout<<" ###Calculating Dos....###" <<endl;
461464 cout<<" =========================" <<endl;
@@ -479,6 +482,7 @@ void ESolver_SDFT_PW:: caldos( const int nche_dos, const double sigmain, const d
479482 int nchip_new = ceil ((double )this ->stowf .nchip_max / npart);
480483 allorderchi = new std::complex <double > [nchip_new * npwx * nche_dos];
481484 }
485+ ModuleBase::timer::tick (this ->classname ," Tracepoly" );
482486 cout<<" 1. TracepolyA:" <<endl;
483487 for (int ik = 0 ;ik < nk;ik++)
484488 {
@@ -532,32 +536,32 @@ void ESolver_SDFT_PW:: caldos( const int nche_dos, const double sigmain, const d
532536 }
533537 if (stoiter.method == 2 ) delete[] allorderchi;
534538
535- string dosfile = GlobalV::global_out_dir+" DOS1_smearing.dat" ;
536- ofstream ofsdos (dosfile.c_str ());
539+ ofstream ofsdos;
537540 int ndos = int ((emax-emin) / de)+1 ;
538- double *dos = new double [ndos];
539- ModuleBase::GlobalFunc::ZEROS (dos,ndos);
540541 stoiter.stofunc .sigma = sigmain / ModuleBase::Ry_to_eV;
541- double sum = 0 ;
542- double maxerror = 0 ;
543- ofsdos<<setw (8 )<<" ## E(eV) " <<setw (20 )<<" dos(eV^-1)" <<setw (20 )<<" sum" <<setw (20 )<<" Error(eV^-1)" <<endl;
542+ ModuleBase::timer::tick (this ->classname ," Tracepoly" );
543+
544544 cout<<" 2. Dos:" <<endl;
545+ ModuleBase::timer::tick (this ->classname ," DOS Loop" );
545546 int n10 = ndos/10 ;
546547 int percent = 10 ;
548+ double *sto_dos = new double [ndos];
549+ double *ks_dos = new double [ndos];
550+ double *error = new double [ndos];
547551 for (int ie = 0 ; ie < ndos; ++ie)
548552 {
549- double KS_dos = 0 ;
550- double sto_dos = 0 ;
553+ double tmpks = 0 ;
554+ double tmpsto = 0 ;
551555 stoiter.stofunc .targ_e = (emin + ie * de) / ModuleBase::Ry_to_eV;
552556 if (stoiter.method == 1 )
553557 {
554558 che.calcoef_real (&stoiter.stofunc , &Sto_Func<double >::ngauss);
555- sto_dos = BlasConnector::dot (nche_dos,che.coef_real ,1 ,spolyv,1 );
559+ tmpsto = BlasConnector::dot (nche_dos,che.coef_real ,1 ,spolyv,1 );
556560 }
557561 else
558562 {
559563 che.calcoef_real (&stoiter.stofunc , &Sto_Func<double >::nroot_gauss);
560- sto_dos = stoiter.vTMv (che.coef_real ,spolyv,nche_dos);
564+ tmpsto = stoiter.vTMv (che.coef_real ,spolyv,nche_dos);
561565 }
562566 if (GlobalV::NBANDS > 0 )
563567 {
@@ -566,46 +570,66 @@ void ESolver_SDFT_PW:: caldos( const int nche_dos, const double sigmain, const d
566570 double *en=&(this ->pelec ->ekb (ik, 0 ));
567571 for (int ib = 0 ; ib < GlobalV::NBANDS; ++ib)
568572 {
569- KS_dos += stoiter.stofunc .gauss (en[ib]) * GlobalC::kv.wk [ik] / 2 ;
573+ tmpks += stoiter.stofunc .gauss (en[ib]) * GlobalC::kv.wk [ik] / 2 ;
570574 }
571575 }
572576 }
573- KS_dos /= GlobalV::NPROC_IN_POOL;
574- #ifdef __MPI
575- MPI_Allreduce (MPI_IN_PLACE, &KS_dos, 1 , MPI_DOUBLE, MPI_SUM , STO_WORLD);
576- MPI_Allreduce (MPI_IN_PLACE, &sto_dos, 1 , MPI_DOUBLE, MPI_SUM , MPI_COMM_WORLD);
577- #endif
578- double tmpre = 0 ;
577+ tmpks /= GlobalV::NPROC_IN_POOL;
578+
579+ double tmperror = 0 ;
579580 if (stoiter.method == 1 )
580581 {
581- tmpre = che.coef_real [nche_dos-1 ] * spolyv[nche_dos-1 ];
582+ tmperror = che.coef_real [nche_dos-1 ] * spolyv[nche_dos-1 ];
582583 }
583584 else
584585 {
585586 const int norder = nche_dos;
586587 double last_coef = che.coef_real [norder-1 ];
587588 double last_spolyv = spolyv[norder*norder - 1 ];
588- tmpre = last_coef *(BlasConnector::dot (norder,che.coef_real ,1 ,spolyv+norder*(norder-1 ),1 )
589+ tmperror = last_coef *(BlasConnector::dot (norder,che.coef_real ,1 ,spolyv+norder*(norder-1 ),1 )
589590 + BlasConnector::dot (norder,che.coef_real ,1 ,spolyv+norder-1 ,norder)-last_coef*last_spolyv);
590591 }
591- #ifdef __MPI
592- MPI_Allreduce (MPI_IN_PLACE, &tmpre, 1 , MPI_DOUBLE, MPI_SUM , MPI_COMM_WORLD);
593- #endif
594- if (maxerror < tmpre) maxerror = tmpre;
595- dos[ie] = (KS_dos + sto_dos) / ModuleBase::Ry_to_eV;
596- sum += dos[ie];
597- ofsdos <<setw (8 )<< emin + ie * de <<setw (20 )<<dos[ie]<<setw (20 )<<sum * de <<setw (20 ) <<tmpre <<endl;
592+
598593 if (ie%n10 == n10 -1 )
599594 {
600595 cout<<percent<<" %" <<" " ;
601596 percent+=10 ;
602597 }
598+ sto_dos[ie] = tmpsto;
599+ ks_dos[ie] = tmpks;
600+ error[ie] = tmperror;
603601 }
604- cout<<endl;
605- cout<<" Finish DOS" <<endl;
606- cout<<scientific<<" DOS max Chebyshev Error: " <<maxerror<<endl;
607- delete[] dos;
602+ #ifdef __MPI
603+ MPI_Allreduce (MPI_IN_PLACE, ks_dos, ndos, MPI_DOUBLE, MPI_SUM , STO_WORLD);
604+ MPI_Allreduce (MPI_IN_PLACE, sto_dos, ndos, MPI_DOUBLE, MPI_SUM , MPI_COMM_WORLD);
605+ MPI_Allreduce (MPI_IN_PLACE, error, ndos, MPI_DOUBLE, MPI_SUM , MPI_COMM_WORLD);
606+ #endif
607+ if (GlobalV::MY_RANK == 0 )
608+ {
609+ string dosfile = GlobalV::global_out_dir+" DOS1_smearing.dat" ;
610+ ofsdos.open (dosfile.c_str ());
611+ double maxerror = 0 ;
612+ double sum = 0 ;
613+ ofsdos<<setw (8 )<<" ## E(eV) " <<setw (20 )<<" dos(eV^-1)" <<setw (20 )<<" sum" <<setw (20 )<<" Error(eV^-1)" <<endl;
614+ for (int ie = 0 ; ie < ndos ; ++ie)
615+ {
616+ double tmperror = abs (error[ie]);
617+ if (maxerror < tmperror) maxerror = tmperror;
618+ double dos = (ks_dos[ie] + sto_dos[ie]) / ModuleBase::Ry_to_eV;
619+ sum += dos;
620+ ofsdos <<setw (8 )<< emin + ie * de <<setw (20 )<< dos <<setw (20 )<< sum * de <<setw (20 )<< tmperror <<endl;
621+ }
622+ cout<<endl;
623+ cout<<" Finish DOS" <<endl;
624+ cout<<scientific<<" DOS max Chebyshev Error: " <<maxerror<<endl;
625+ ofsdos.close ();
626+ }
627+ delete[] sto_dos;
628+ delete[] ks_dos;
629+ delete[] error;
608630 delete[] spolyv;
631+ ModuleBase::timer::tick (this ->classname ," DOS Loop" );
632+ ModuleBase::timer::tick (this ->classname ," caldos" );
609633 return ;
610634}
611635
0 commit comments