block_matrix_t atomic_density_matrix(atom_diag const& atom, double beta) {
 double z = partition_function(atom, beta);
 int n_blocks = atom.n_blocks();
 block_matrix_t dm(n_blocks);
 for (int bl = 0; bl < n_blocks; ++bl) {
  int bl_size = atom.get_block_dim(bl);
  dm[bl] = matrix<h_scalar_t>(bl_size, bl_size);
  for (int u = 0; u < bl_size; ++u) {
   for (int v = 0; v < bl_size; ++v) {
    dm[bl](u, v) = (u == v) ? std::exp(-beta * atom.get_eigenvalue(bl, u)) / z : 0;
   }
  }
 }
 return dm;
}
block_gf<imtime> atomic_gf(atom_diag const& atom, double beta, std::map<std::string, indices_t> const& gf_struct,
                                      int n_tau, std::vector<std::pair<int, int>> const& excluded_states) {

 double z = partition_function(atom, beta);

 std::vector<std::string> block_names;
 std::vector<gf<imtime>> gf_blocks;
 auto const & fops = atom.get_fops();

 auto is_excluded = [&excluded_states](int A, int ia) {
  return std::find(excluded_states.begin(), excluded_states.end(), std::make_pair(A, ia)) != excluded_states.end();
 };

 for (auto const& block : gf_struct) {
  block_names.push_back(block.first);
  int bl_size = block.second.size();
  auto g = gf<imtime>{{beta, Fermion, n_tau}, {bl_size, bl_size}};

  for (int inner_index1 = 0; inner_index1 < bl_size; ++inner_index1)
   for (int inner_index2 = 0; inner_index2 < bl_size; ++inner_index2) {
    int n1 = fops[{block.first, block.second[inner_index1]}]; // linear_index of c
    int n2 = fops[{block.first, block.second[inner_index2]}]; // linear_index of c_dag

    for (int A = 0; A < atom.n_blocks(); ++A) {              // index of the A block. sum over all
     int B = atom.cdag_connection(n2, A);                // index of the block connected to A by operator c_n
     if (B == -1) continue;                             // no matrix element
     if (atom.c_connection(n1, B) != A) continue; //
     for (int ia = 0; ia < atom.get_block_dim(A); ++ia)
      for (int ib = 0; ib < atom.get_block_dim(B); ++ib) {
       auto Ea = atom.get_eigenvalue(A, ia);
       auto Eb = atom.get_eigenvalue(B, ib);
       if (is_excluded(A, ia) || is_excluded(B, ib)) continue;
       for (auto tau : g.mesh())
        g[tau](inner_index1, inner_index2) +=
            -atom.cdag_matrix(n2, A)(ib, ia) * atom.c_matrix(n1, B)(ia, ib) * std::exp(-(Eb - Ea) * tau - beta * Ea) / z;
      }
    }
   }
  g.singularity()(1) = 1.0;
  gf_blocks.push_back(std::move(g));
 }

 return make_block_gf(block_names, gf_blocks);
}
Exemple #3
0
/* ==== */
static void
wl_montecarlo(char *struc)
{
  short *pt=NULL;
  move_str m;
  int e,enew,emove,eval_me,status,debug=1;
  long int crosscheck=1000000; /* used for convergence checks */
  long int crosscheck_limit = 100000000000000000;
  double g_b1,g_b2,prob,lnf = 1.;  /* log modification parameter f */
  size_t b1,b2;                    /* indices in g/h corresponding to
				      old/new energies */
  gsl_histogram *gcp=NULL; /* clone of g used during crosscheck output */ 

  eval_me = 1; /* paranoid checking of neighbors against RNAeval */
  if (wanglandau_opt.verbose){
    printf("[[wl_montecarlo()]]\n");
  }
  pt = vrna_pt_get(struc);
  //mtw_dump_pt(pt);
  //char *str = vrna_pt_to_db(pt);
  //printf(">%s<\n",str);
  e = vrna_eval_structure_pt(wanglandau_opt.sequence,pt,P);
  
  /* determine bin where the start structure goes */
  status = gsl_histogram_find(g,(float)e/100,&b1);
  if (status) {
    if (status == GSL_EDOM){
      printf ("error: %s\n", gsl_strerror (status));
    }
    else {fprintf(stderr, "GSL error: gsl_errno=%d\n",status);}
    exit(EXIT_FAILURE);
  }
  printf("%s\n", wanglandau_opt.sequence);
  print_str(stderr,pt);
  printf(" (%6.2f) bin:%d\n",(float)e/100,b1);
  if (wanglandau_opt.verbose){
    fprintf(stderr,"\nStarting MC loop ...\n");
  }
  while (lnf > wanglandau_opt.ffinal) {
    if(wanglandau_opt.debug){
      fprintf(stderr,"\n==================\n");
      fprintf(stderr,"in while: lnf=%8.6f\n",lnf);
      fprintf(stderr,"steps: %d\n",steps);
      fprintf(stderr,"current histogram g:\n");
      gsl_histogram_fprintf(stderr,g,"%6.2f","%30.6f");
      fprintf(stderr,"\n");
      print_str(stderr,pt);
      fprintf(stderr, " (%6.2f) bin:%d\n",(float)e/100,b1);
      /*  mtw_dump_pt(pt); */
    }
    /* make a random move */
    m = get_random_move_pt(wanglandau_opt.sequence,pt);
    /* compute energy difference for this move */
    emove = vrna_eval_move_pt(pt,s0,s1,m.left,m.right,P);
    /* evaluate energy of the new structure */
    enew = e + emove;
    if(wanglandau_opt.debug){
      fprintf(stderr,
	      "random move: left %i right %i enew(%6.4f)=e(%6.4f)+emove(%6.4f)\n",
	      m.left,m.right,(float)enew/100,(float)e/100,(float)emove/100);
    }

    /* ensure the new energy is within sampling range */
    if ((float)enew/100 >= wanglandau_opt.max){
      fprintf(stderr,
	      "New structure has energy %6.2f >= %6.2f (upper energy bound)\n",
	      (float)enew/100,wanglandau_opt.max);
      fprintf(stderr,"Please increase --bins or adjust --max! Exiting ...\n");
      exit(EXIT_FAILURE);
    }
    /* determine bin where the new structure goes */
    status = gsl_histogram_find(g,(float)enew/100,&b2);
    if (status) {
      if (status == GSL_EDOM){
	printf ("error: %s\n", gsl_strerror (status));
      }
      else {fprintf(stderr, "GSL error: gsl_errno=%d\n",status);}
      exit(EXIT_FAILURE);
    }

    steps++;  /* # of MC steps performed so far */

    /* lookup current values for bins b1 and b2 */
    g_b1 = gsl_histogram_get(g,b1);
    g_b2 = gsl_histogram_get(g,b2);
      
    /* core MC steps */
    prob = MIN2(exp(g_b1 - g_b2), 1.0);
    rnum =  gsl_rng_uniform (r);
    
    if ((prob == 1 || (rnum <= prob)) ) { /* accept & apply the move */
      apply_move_pt(pt,m);
      if(wanglandau_opt.debug){
	print_str(stderr,pt);
	fprintf(stderr, " %6.2f bin:%d [A]\n", (float)enew/100,b2);
      }
      b1 = b2;
      e = enew;
    }
    else { /* reject the move */
      if(wanglandau_opt.debug){
	print_str(stderr,pt);
	fprintf(stderr, " (%6.2f) bin:%d [R]\n", (float)enew/100,b2);
       }
    }
    
    /* update histograms g and h */
    if(wanglandau_opt.truedosbins_given && b2 <= wanglandau_opt.truedosbins){
      /* do not update if b2 <= truedosbins, i.e. keep true DOS values
	 in those bins */
      if (wanglandau_opt.debug){
	fprintf(stderr, "NOT UPDATING bin %d\n",b1);
      }
    } else{
      if(wanglandau_opt.debug){
	fprintf(stderr, "UPDATING bin %d\n",b1); 
      }
      status = gsl_histogram_increment(h,(float)e/100);
      status = gsl_histogram_accumulate(g,(float)e/100,lnf);
    }
    maxbin = MAX2(maxbin,(int)b1);
   
    // stuff that can be skipped 
    /*
      printf ("performed move l:%4d r:%4d\t Energy +/- %6.2f\n",m.left,m.right,(float)emove/100);
      print_str(stderr,pt);printf(" %6.2f bin:%d\n",(float)enew/100,b2);
      e = vrna_eval_structure_pt(wanglandau_opt.sequence,pt,P);
      if (eval_me == 1 && e != enew){
      fprintf(stderr, "energy evaluation against vrna_eval_structure_pt() mismatch... HAVE %6.2f != %6.2f (SHOULD BE)\n",(float)enew/100, (float)e/100);
      exit(EXIT_FAILURE);
      }
      print_str(stderr,pt);printf(" %6.2f\n",(float)e/100);
    */
    // end of stuff that can be skipped

    /* output DoS every x*10^(1/4) steps, starting with x=10^6 (we
       used this fopr comparing perfomance and convergence of
       different DoS sampling methods */
    if((steps % crosscheck == 0) && (crosscheck <= crosscheck_limit)){
      fprintf(stderr,"# crosscheck reached %li steps ",crosscheck);
      gcp = gsl_histogram_clone(g);
      if(wanglandau_opt.verbose){
	fprintf(stderr,"## gcp before scaling\n");
	gsl_histogram_fprintf(stderr,gcp,"%6.2f","%30.6f");
      }
      scale_dos(gcp); /* scale estimated g; make ln(g[0])=0 */
      if(wanglandau_opt.verbose){
	fprintf(stderr,"## gcp after scaling\n");
	gsl_histogram_fprintf(stderr,gcp,"%6.2f","%30.6f");
      }
      double Z = partition_function(gcp);
      output_dos(gcp,'s');
      fprintf(stderr, "Z=%10.4g\n", Z);
      crosscheck *= (pow(10, 1.0/4.0));
      gsl_histogram_free(gcp);
      fprintf(stderr,"->  new crosscheck will be performed at %li steps\n", crosscheck);
    }
    
    if(steps % wanglandau_opt.checksteps == 0) {
      if( histogram_is_flat(h) ) {
	lnf /= 2;
	fprintf(stderr,"# steps=%20li | f=%12g | histogram is FLAT\n",
		steps,lnf);
	gsl_histogram_reset(h);
      }
      else {
	fprintf(stderr, "# steps=%20li | f=%12g | histogram is NOT FLAT\n",
		steps,lnf);
      }
      output_dos(g,'l');
    }
    
    /* stop criterion */
    if(steps % wanglandau_opt.steplimit == 0){
      fprintf(stderr,"maximun number of MC steps (%li) reached, exiting ...",
	      wanglandau_opt.steplimit);
      break;
    }

  } /* end while */
  free(pt); 
  return;
}