def good2n(nmax, coefs_list, ref_coefs, threshold=0.90, outfile=''):

    ## cc=0.90 is equivalent to 5% mutation in real space at nmax<=10
    max_indx = math.nlm_array(nmax).nlm().size()
    for nn in range(nmax, 1, -1):
        min_indx = math.nlm_array(nn - 1).nlm().size()
        #coef_0 = ref_coefs[min_indx:max_indx]
        coef_0 = ref_coefs[0:max_indx]
        mean_0 = abs(ref_coefs[0])
        sigma_0 = flex.sum(flex.norm(coef_0)) - mean_0**2
        sigma_0 = smath.sqrt(sigma_0)
        cc_array = flex.double()
        #out = open(outfile,"w")
        for coef in coefs_list:
            #coef_1 = coef[min_indx:max_indx]
            coef_1 = coef[0:max_indx]
            mean_1 = abs(coef[0])
            sigma_1 = flex.sum(flex.norm(coef_1)) - mean_1**2
            sigma_1 = smath.sqrt(sigma_1)
            cov_01 = abs(flex.sum(coef_0 * flex.conj(coef_1)))
            cov_01 = cov_01 - mean_0 * mean_1
            this_cc = cov_01 / sigma_1 / sigma_0
            cc_array.append(this_cc)
            out = open(outfile, "a")
            print >> out, this_cc
            out.close()
            print this_cc
        mean_cc = flex.mean(cc_array)
        out = open(outfile, "a")
        print >> out, "level n: ", nn, mean_cc
        out.close()
        print "level n: ", nn, mean_cc
        if (mean_cc >= threshold): return nn
        max_indx = min_indx
    return nn
Exemplo n.º 2
0
def test_direct_sum():

    # correct values
    p = pdb.input(source_info='string', lines=test_pdb)
    x = p.xray_structure_simple()
    for s in x.scatterers():
        s.set_use_u(False, False)

    fc = x.structure_factors(anomalous_flag=False,
                             d_min=2.0,
                             algorithm='direct').f_calc()
    fcd = fc.data()
    indices = fc.indices()

    # test values
    xyz = x.sites_frac()
    h = flex.vec3_double(len(indices))
    fm = matrix.sqr(
        p.crystal_symmetry().unit_cell().fractionalization_matrix())
    om = matrix.sqr(
        p.crystal_symmetry().unit_cell().orthogonalization_matrix())
    for i in xrange(len(indices)):
        h[i] = fm * indices[i]
    sr = x.scattering_type_registry()
    st = x.scattering_types()

    sg = p.crystal_symmetry().space_group()
    r = flex.double()
    t = flex.vec3_double(len(sg))
    for i in xrange(len(sg)):
        r_i = om * matrix.sqr(sg[i].r().as_double())
        for a in r_i:
            r.append(a)
        t[i] = om * matrix.col(sg[i].t().as_double())

    bls = flex.double(len(xyz), 0.0)
    amplitudes = cuda_direct_summation()
    amplitudes.add(st, xyz, bls, h, r, t, sr)
    amplitudes = amplitudes.get_sum()

    cpu_i = flex.norm(fcd)
    gpu_i = flex.norm(amplitudes)

    mean = 0.0
    for i in xrange(len(cpu_i)):
        e = math.fabs(cpu_i[i] - gpu_i[i]) / cpu_i[i]
        mean += e
    mean = mean / (len(cpu_i))
    assert (mean < 1.0e-3)
    def build_image(self, n_photons=1e20, coherent=True):
        """
    From the International Tables of Crystallography (2006). Vol. C, Chapter 2.6,
    the scattered intensity from one electron is given the by Thomson formula

                    r_e^2   / 1 + cos^2(2 theta) \
         I_e = I_0 ------- | -------------------- |          (Eq 2.6.1.1)
                     r^2    \         2          /

    where r_e is the classical radius of the electron (2.818 x 10^(-15) m),
    r is the distance between the sample and the detector, and 2 theta is
    the scattering angle.  The total scattering intensity is thus I_e
    multiplied by the total number of scattering electrons.
    """
        # cache Miller indices
        if (self.cached_h is None):
            self.cached_h = self.image_base.get_corner_h()
            self.cached_q = self.image_base.get_center_q()

        # calculate intensities
        if (gpu_available):
            if (coherent):
                self.sum_structure_factors_gpu()
                self.intensities = flex.norm(self.structure_factors)
            else:
                self.incoherent_sum_structure_factors_gpu()
        else:
            if (coherent):
                self.sum_structure_factors_cpu()
                self.intensities = flex.norm(self.structure_factors)
            else:
                self.incoherent_sum_structure_factors_cpu()

        # scale intensities
        r_e = 2.818e-5  # radius of electron in Angstroms
        d = self.image_base.get_ewald_sphere().get_distance(
        )  # detector distance
        i000 = 0.0
        for i in xrange(len(self.structures)):
            i000 += self.structures.species[i].scattering_type_registry.\
                    sum_of_scattering_factors_at_diffraction_angle_0() *\
                    self.structures.species[i].n_copies
        if (approx_equal(i000, 0.0, out=None)):
            i000 = 1.0
        scale = self.structures.total_electrons / i000 * n_photons * (
            r_e * r_e) / (d * d)
        self.intensities = scale * self.intensities
        return self.intensities
    def build_image(self, n_photons=1e20, coherent=True):
        """
    From the International Tables of Crystallography (2006). Vol. C, Chapter 2.6,
    the scattered intensity from one electron is given the by Thomson formula

                    r_e^2   / 1 + cos^2(2 theta) \
         I_e = I_0 ------- | -------------------- |          (Eq 2.6.1.1)
                     r^2    \         2          /

    where r_e is the classical radius of the electron (2.818 x 10^(-15) m),
    r is the distance between the sample and the detector, and 2 theta is
    the scattering angle.  The total scattering intensity is thus I_e
    multiplied by the total number of scattering electrons.
    """
        r_e = 2.818e-5  # radius of electron in Angstroms
        d = self.image_composer.get_distance()  # detector distance
        if (self.structure_factors is None):
            self.sum_structure_factors_gpu()
        if (coherent):
            intensities = flex.norm(self.structure_factors)
        else:
            intensities = self.structure_factors
        bc = self.image_composer.get_beam_center()
        ds = self.image_composer.get_detector_size()
        i000 = intensities[bc[1] * ds[0] + bc[0]]
        if (approx_equal(i000, 0.0, out=None)):
            i000 = 1.0
        scale = self.structures.total_electrons / i000 * n_photons * (
            r_e * r_e) / (d * d)
        intensities = scale * intensities
        image = self.image_composer.build_image(intensities)
        return image
    def incoherent_sum_structure_factors_cpu(self):

        self.intensities = flex.complex_double(len(self.cached_h),
                                               complex(0.0, 0.0))

        # sum structure factors incoherently
        for i_structure in xrange(len(self.structures.species)):
            for j_orientation in xrange(
                    self.structures.species[i_structure].n_copies):
                tmp_xyz = self.structures.species[i_structure].xyz.deep_copy()
                translation = flex.double(
                    self.structures.translations[i_structure][j_orientation])
                rotation = matrix.sqr(
                    self.structures.rotations[i_structure][j_orientation])
                for i in xrange(len(tmp_xyz)):
                    xyz = flex.double(rotation * tmp_xyz[i])
                    tmp_xyz[i] = tuple(xyz + translation)
                structure_factors = direct_sum_structure_factors\
                                    (self.structures.species[i_structure].
                                     scattering_types,
                                     tmp_xyz,
                                     self.structures.species[i_structure].
                                     boundary_layer_scaling_factors,
                                     self.cached_h,
                                     self.structures.species[i_structure].
                                      scattering_type_registry)
                self.intensities += flex.norm(structure_factors)
Exemplo n.º 6
0
  def pair_align(self):
    ms = flex.double()
    ss = flex.double()
    tmp_nlm_array = math.nlm_array( self.nmax )
    for coef in self.finals:
      mean = abs( coef[0] )
      var = flex.sum( flex.norm( coef ) )
      sigma = smath.sqrt( var - mean*mean )
      ms.append( mean )
      ss.append( sigma)

    grids = flex.grid(self.n_trial, self.n_trial)
    self.cc_array=flex.double( grids, 1.0 )
    for ii in range( self.n_trial ):
      self.nlm_array.load_coefs( self.nlm, self.finals[ii] )
      for jj in range( ii ):
        tmp_nlm_array.load_coefs( self.nlm, self.finals[jj] )
        cc = fft_align.align( self.nlm_array, tmp_nlm_array, nmax=self.nmax, refine=True ).best_score
        cc = (cc-ms[ii]*ms[jj])/(ss[ii]*ss[jj])
        self.cc_array[(ii,jj)]=cc
        self.cc_array[(jj,ii)]=cc

    outfile = self.prefix+"pair.cc"
    comment = "# electron density correlation coefficient, < rho_1(r)*rho_2(r) >"
    out=open(outfile, 'w')
    print>>out, comment
    for ii in range(1,self.n_trial+1):
      print>>out,"%6d"%ii,
    print>>out, "   average"

    for ii in range(self.n_trial):
      for jj in range(self.n_trial):
        print>>out,"%6.3f"%self.cc_array[(ii,jj)],
      print>>out, flex.mean( self.cc_array[ii*self.n_trial:(ii+1)*self.n_trial] )
    out.close()
def test_direct_summation():

  # correct values
  p = pdb.input(source_info='string',lines=test_pdb)
  x = p.xray_structure_simple()
  for s in x.scatterers():
    s.set_use_u(False,False)

  fc = x.structure_factors(anomalous_flag=False,d_min=2.0,
                           algorithm='direct').f_calc()
  fcd = fc.data()
  indices = fc.indices()

  # test values
  xyz = x.sites_frac()
  h = flex.vec3_double(len(indices))
  fm = matrix.sqr(p.crystal_symmetry().unit_cell().fractionalization_matrix())
  om = matrix.sqr(p.crystal_symmetry().unit_cell().orthogonalization_matrix())
  for i in xrange(len(indices)):
    h[i] = fm * indices[i]
  sr = x.scattering_type_registry()
  st = x.scattering_types()

  sg = p.crystal_symmetry().space_group()
  r = flex.double()
  t = flex.vec3_double(len(sg))
  for i in xrange(len(sg)):
    r_i = om * matrix.sqr(sg[i].r().as_double())
    for a in r_i:
      r.append(a)
    t[i] = om * matrix.col(sg[i].t().as_double())

  bls = flex.double(len(xyz),0.0)
  amplitudes = direct_summation()
  amplitudes.add(st,xyz,bls,h,r,t,sr,False)
  amplitudes = amplitudes.get_sum()

  cpu_i = flex.norm(fcd)
  gpu_i = flex.norm(amplitudes)

  mean = 0.0
  for i in xrange(len(cpu_i)):
    e = math.fabs(cpu_i[i] - gpu_i[i])/cpu_i[i]
    mean += e
  mean = mean/(len(cpu_i))
  assert(mean < 1.0e-3)
Exemplo n.º 8
0
def get_mean_sigma(nlm_array):
    t1 = time.time()
    coef = nlm_array.coefs()
    mean = abs(coef[0])
    var = flex.sum(flex.norm(coef))
    sigma = smath.sqrt(var - mean * mean)
    t2 = time.time()
    print("get_mean_sigma time used:", t2 - t1)
    return mean, sigma
Exemplo n.º 9
0
  def refine(self, trial):
    print "--------------Trial %d-----------------"%trial, time.ctime()
    self.working_model = self.start_model.deep_copy()
    self.nlm_coefs = self.start_nlm_coefs.deep_copy()
    self.best_nlm_coefs = self.start_nlm_coefs.deep_copy()
    self.best_blq = self.start_blq.deep_copy()
    self.lowest_score = self.start_score
    init_scores = flex.double()
    while( init_scores.size() < 10 ):
      if(self.modify()):
        init_scores.append( self.target() )
    mean = flex.mean( init_scores )
    self.deltaS = smath.sqrt( flex.sum(flex.pow2(init_scores-mean) )/init_scores.size() )
    self.T = self.deltaS * 20
    self.nsteps = 500
    self.score = mean
    self.working_model = self.start_model.deep_copy()
    while( True ): #self.T > self.deltaS/10.0):
      self.n_reject_this_round = 0
      self.n_accept_this_round = 0
      for ii in range( self.nsteps ):
        self.move()
      self.n_moves_this_round = self.n_reject_this_round + self.n_accept_this_round
      print "Number of moves: %d(out of %d)"%(self.n_moves_this_round, self.nsteps)
      print "Number of Accept/Reject: %d/%d"%(self.n_accept_this_round, self.n_reject_this_round)
      if( self.n_reject_this_round >= self.n_moves_this_round*0.9 ):
        print "Too Many rejections (%d), quit at temperature (%f)"%(self.n_reject_this_round, self.T)
        break
      self.T = self.T*0.9

    self.nlm_array.load_coefs( self.nlm, self.best_nlm_coefs )
    best_blq = self.blq_calculator.get_all_blq( self.nlm_array )
    best_blq = best_blq/best_blq[0]
    out = open(self.prefix+str(trial)+'_final.blq', 'w')
    self.data.print_out(data=best_blq,out=out)
    out.close()
    print "total number of moves %d"%self.counter
    print "total number of accepted moves %d"%self.n_accept

    if (self.pdb_nlm is not None):
      align_obj = fft_align.align(self.pdb_nlm, self.nlm_array, nmax=self.nmax, refine=True)
      mean = abs( self.best_nlm_coefs[0] )
      var = flex.sum( flex.norm( self.best_nlm_coefs ) )
      sigma = smath.sqrt( var - mean*mean )
      cc = align_obj.best_score
      cc = ( cc - mean*self.pdb_m ) / ( sigma*self.pdb_s )
      print "C.C. (PDB, trial%6d) = %8.5f, Score = %8.5f"%(trial, cc, self.lowest_score)
      self.best_nlm_coefs = align_obj.moving_nlm.coefs()

    reconst_model = self.reconst_model( self.best_nlm_coefs )
    xplor_map_type( reconst_model, self.np_on_grid, self.rmax, file_name=self.prefix+str(trial)+'_final_rbt.xplor')
    xplor_map_type( self.best_model, self.np_on_grid, self.rmax, file_name=self.prefix+str(trial)+'_final.xplor')
    print "-----------End of Trial %d--------------"%trial, time.ctime()
Exemplo n.º 10
0
  def refine(self, trial):
    print "--------------Trial %d-----------------"%trial, time.ctime()
    self.working_model = self.start_model.deep_copy()
    self.nlm_coefs = self.start_nlm_coefs.deep_copy()
    self.best_nlm_coefs = self.start_nlm_coefs.deep_copy()
    self.best_i = self.start_i.deep_copy()
    self.lowest_score = self.start_score
    init_scores = flex.double()
    for ii in range(10):
      self.modify()
      init_scores.append( self.target() )
    mean = flex.mean( init_scores )
    self.deltaS = smath.sqrt( flex.sum(flex.pow2(init_scores-mean) )/10.0 )
    self.T = self.deltaS * 100
    self.nsteps = 200
    self.score = mean
    self.working_model = self.start_model.deep_copy()
    while( self.T > self.deltaS/2.0):
      self.n_reject = 0
      for ii in range( self.nsteps ):
        self.move()
      print "Number of Accept/Reject: %d/%d"%(self.nsteps-self.n_reject, self.n_reject)
      if( self.n_reject > self.nsteps*0.9 ):
        print "Too Many rejections (%d), quit at temperature (%f)"%(self.n_reject, self.T)
        break
      self.T = self.T*0.9

    out = open(self.prefix+str(trial)+'_final.iq', 'w')
    self.nlm_array.load_coefs( self.nlm, self.best_nlm_coefs )
    best_i = self.zm.calc_intensity_nlm( self.nlm_array )
    best_i = best_i/best_i[0]*self.scale_2_expt
    for qq,ic,io in zip( self.data.q, best_i, self.data.i*self.scale_2_expt):
      print>>out, qq, ic, io
    out.close()
    print "total number of moves %d"%self.counter
    print "total number of accepted moves %d"%self.n_accept

    if (self.pdb_nlm is not None):
      align_obj = fft_align.align(self.pdb_nlm, self.nlm_array, nmax=self.nmax, refine=True)
      mean = abs( self.best_nlm_coefs[0] )
      var = flex.sum( flex.norm( self.best_nlm_coefs ) )
      sigma = smath.sqrt( var - mean*mean )
      cc = align_obj.best_score
      cc = ( cc - mean*self.pdb_m ) / ( sigma*self.pdb_s )
      print "C.C. (PDB, trial%6d) = %8.5f, Score = %8.5f"%(trial, cc, self.lowest_score)
      self.best_nlm_coefs = align_obj.moving_nlm.coefs()

    reconst_model = self.reconst_model( self.best_nlm_coefs )
    xplor_map_type( reconst_model, self.np_on_grid, self.rmax, file_name=self.prefix+str(trial)+'_final_rbt.xplor')
    xplor_map_type( self.best_model, self.np_on_grid, self.rmax, file_name=self.prefix+str(trial)+'_final.xplor')
    print "-----------End of Trial %d--------------"%trial, time.ctime()
    def incoherent_sum_structure_factors_gpu(self):

        self.intensities = flex.double(len(self.cached_h), 0.0)

        # sum structure factors incoherently
        for i in xrange(len(self.structures)):
            for rt in xrange(len(self.structures.rotations[i])):
                sf = cuda_direct_summation()
                rot = flex.double()
                trans = flex.vec3_double()
                for j in xrange(len(self.structures.rotations[i][rt])):
                    rot.append(self.structures.rotations[i][rt][j])
                trans.append(self.structures.translations[i][rt])
                sf.add(self.structures.species[i].scattering_types,
                       self.structures.species[i].xyz,
                       self.structures.species[i].\
                       boundary_layer_scaling_factors,
                       self.cached_h,rot,trans,
                       self.structures.species[i].\
                       scattering_type_registry)
                self.intensities += flex.norm(sf.get_sum())
Exemplo n.º 12
0
def get_mean_sigma( nlm_array ):
  coef = nlm_array.coefs()
  mean = abs( coef[0] )
  var = flex.sum( flex.norm(coef) )
  sigma = math.sqrt( var-mean*mean )
  return mean, sigma
Exemplo n.º 13
0
  def scan( self ):
    fft = fftpack.complex_to_complex_2d( self.ngrid, self.ngrid )
    inversion = False
    for beta in self.beta:
      self.cc_obj.set_beta( beta )
      mm = self.cc_obj.mm_coef(0,inversion)
      if( self.pad > 0):
        mm = self.cc_obj.mm_coef(self.pad, inversion)
      fft_input= mm
      scores = fft.backward( fft_input ).as_1d()
      self.scores = self.scores.concatenate( -flex.norm( scores )  )
    self.best_indx = flex.min_index( self.scores )
    self.best_score = math.sqrt( -self.scores[ self.best_indx ])


    if self.check_inversion:
    ### Inversion of the Spherical Harmonics ###
      inversion = True
      inversion_scores = flex.double()
      for beta in self.beta:
        self.cc_obj.set_beta( beta )
        mm = self.cc_obj.mm_coef(0,inversion)
        if( self.pad > 0):
          mm = self.cc_obj.mm_coef(self.pad, inversion)
        fft_input= mm.deep_copy()
        scores = fft.backward( fft_input ).as_1d()
        inversion_scores = inversion_scores.concatenate( -flex.norm( scores )  )
      inv_best_indx = flex.min_index( inversion_scores )
      inv_best_score = math.sqrt(-inversion_scores[ inv_best_indx ] )

      if( inv_best_score < self.best_score ):
        self.score = inversion_scores
        self.best_indx = inv_best_indx
        self.best_score = inv_best_score
        self.inversion =  True
      else:
        self.inversion = False



    b=self.best_indx//(self.ngrid*self.ngrid)
    a=(self.best_indx - self.ngrid*self.ngrid*b ) // self.ngrid
    g=self.best_indx - self.ngrid*self.ngrid*b - self.ngrid*a

    b = self.beta[b]
    g = math.pi*2.0 *( float(g)/(self.ngrid-1) )
    a = math.pi*2.0 *( float(a)/(self.ngrid-1) )

    self.best_ea = (a, b, g )

    self.find_top( self.topn )
    if( self.refine ):
      self.refined = []
      self.refined_moving_nlm = []
      self.refined_score = flex.double()
      for t in self.top_align:
        r = self.run_simplex( t )
        self.refined.append ( r )
        self.refined_score.append( self.get_cc( self.target( r ) ) )
        self.refined_moving_nlm.append(  self.cc_obj.rotate_moving_obj( r[0],r[1], r[2], self.inversion )  )
      orders=flex.sort_permutation( self.refined_score, True )
      self.best_score = -self.refined_score[orders[0]]


# show the refined results
      if( self.show_result ):
        print("refined results:")
        for ii in range( self.topn ):
          o = orders[ii]
          o = ii
          print(ii, ":", list( self.refined[o] ), ":", self.refined_score[o])
      ea = self.refined[ orders[0] ]
      self.best_ea = (ea[0], ea[1], ea[2] )
      self.moving_nlm = self.cc_obj.rotate_moving_obj( ea[0],ea[1], ea[2], self.inversion )
Exemplo n.º 14
0
    def scan(self):
        fft = fftpack.complex_to_complex_2d(self.ngrid, self.ngrid)
        inversion = False
        for beta in self.beta:
            self.cc_obj.set_beta(beta)
            mm = self.cc_obj.mm_coef(0, inversion)
            if self.pad > 0:
                mm = self.cc_obj.mm_coef(self.pad, inversion)
            fft_input = mm
            scores = fft.backward(fft_input).as_1d()
            self.scores = self.scores.concatenate(-flex.norm(scores))
        self.best_indx = flex.min_index(self.scores)
        self.best_score = smath.sqrt(-self.scores[self.best_indx])

        if self.check_inversion:
            ### Inversion of the Spherical Harmonics ###
            inversion = True
            inversion_scores = flex.double()
            for beta in self.beta:
                self.cc_obj.set_beta(beta)
                mm = self.cc_obj.mm_coef(0, inversion)
                if self.pad > 0:
                    mm = self.cc_obj.mm_coef(self.pad, inversion)
                fft_input = mm.deep_copy()
                scores = fft.backward(fft_input).as_1d()
                inversion_scores = inversion_scores.concatenate(-flex.norm(scores))
            inv_best_indx = flex.min_index(inversion_scores)
            inv_best_score = smath.sqrt(-inversion_scores[inv_best_indx])

            if inv_best_score < self.best_score:
                self.score = inversion_scores
                self.best_indx = inv_best_indx
                self.best_score = inv_best_score
                self.inversion = True
            else:
                self.inversion = False

        b = self.best_indx // (self.ngrid * self.ngrid)
        a = (self.best_indx - self.ngrid * self.ngrid * b) // self.ngrid
        g = self.best_indx - self.ngrid * self.ngrid * b - self.ngrid * a

        b = self.beta[b]
        g = smath.pi * 2.0 * (float(g) / (self.ngrid - 1))
        a = smath.pi * 2.0 * (float(a) / (self.ngrid - 1))

        self.best_ea = (a, b, g)

        self.find_top(self.topn)
        if self.refine:
            self.refined = []
            self.refined_moving_nlm = []
            self.refined_score = flex.double()
            for t in self.top_align:
                r = self.run_simplex(t)
                self.refined.append(r)
                self.refined_score.append(self.get_cc(self.target(r)))
                self.refined_moving_nlm.append(self.cc_obj.rotate_moving_obj(r[0], r[1], r[2], self.inversion))
            orders = flex.sort_permutation(self.refined_score, True)
            self.best_score = -self.refined_score[orders[0]]

            # show the refined results
            if self.show_result:
                print "refined results:"
                for ii in range(self.topn):
                    o = orders[ii]
                    o = ii
                    print ii, ":", list(self.refined[o]), ":", self.refined_score[o]
            ea = self.refined[orders[0]]
            self.best_ea = (ea[0], ea[1], ea[2])
            self.moving_nlm = self.cc_obj.rotate_moving_obj(ea[0], ea[1], ea[2], self.inversion)
Exemplo n.º 15
0
def get_mean_sigma(nlm_array):
    coef = nlm_array.coefs()
    mean = abs(coef[0])
    var = flex.sum(flex.norm(coef))
    sigma = smath.sqrt(var - mean * mean)
    return mean, sigma
Exemplo n.º 16
0
  def __init__(self, x, spans=None, demean=True, detrend=True):

    # Ensure x is copied as it will be changed in-place
    x = flex.double(x).deep_copy()
    n = len(x)

    if detrend:
      t = flex.size_t_range(n).as_double() + 1 - (n + 1)/2
      inv_sumt2 = 1./t.dot(t)
      x = x - flex.mean(x) - x.dot(t) * t * inv_sumt2
    elif demean:
      x -= flex.mean(x)

    # determine frequencies
    stop = ((n - (n % 2)) // 2) + 1
    self.freq = flex.double([i / n for i in range(1, stop)])

    fft = fftpack.real_to_complex(n)
    n = fft.n_real()
    m = fft.m_real()
    x.extend(flex.double(m-n, 0.))
    xf = fft.forward(x)

    # get abs length of complex and normalise by n to get the raw periodogram
    spec = flex.norm(xf) / n

    if spans is None:
      # if not doing smoothing, the spectrum is just the raw periodogram with
      # the uninteresting DC offset removed
      self.spec = spec[1:]
      return

    # for smoothing replace the DC offset term and extend the rest of the
    # sequence by its reverse conjugate, omitting the Nyquist term if it is
    # present
    spec[0] = spec[1]
    end = fft.n_complex() - (fft.n_real() + 1) % 2
    spec.extend(spec[1:end].reversed())

    try:
      # multiple integer spans
      nspans = len(spans)
      m = int(spans[0]) // 2
      multiple = True
    except TypeError:
      # single integer span
      m = int(spans) // 2
      multiple = False

    # construct smoothing kernel
    k = Kernel('modified.daniell', m)

    if multiple:
      for i in range(1, nspans):
        # construct kernel for convolution
        m1 = int(spans[i]) // 2
        k1 = Kernel('modified.daniell', m1)

        # expand coefficients of k to full kernel and zero pad for smoothing
        x1 = flex.double(k1.m, 0.0)
        x1.extend(k.coef.reversed())
        x1.extend(k.coef[1:])
        x1.extend(flex.double(k1.m, 0.0))

        # convolve kernels
        coef = kernapply(x1, k1, circular=True)
        m = len(coef)//2
        coef = coef[m:(len(coef))]
        k = Kernel(coef=coef)

    # apply smoothing kernel
    spec = kernapply(spec, k, circular=True)
    self.spec = spec[1:fft.n_complex()]

    return
def get_intensity_structure(base,FE1_model,FE2_model):
  """Now used stuff learned from test_fmodel_stuff to calculate new data structure for use
     in the program.
  """
  # generate the list of all HKL to be used throughout.
  C2_structures = get_C2_structures()

  energy = 7070.0
  GF_whole7070 = gen_fmodel_with_complex.from_structure(C2_structures[0],energy
                 ).from_parameters(algorithm="fft")
  GF_whole7070 = GF_whole7070.as_P1_primitive()
  f_container7070 = GF_whole7070.get_fmodel()
  Fmodel_whole7070 = f_container7070.f_model
  Fmodel_indices = Fmodel_whole7070.indices() # common structure defines the indices
  MS = Fmodel_whole7070.set() # avoid having to repeatedly calculate indices
  CS = Fmodel_whole7070.crystal_symmetry() # same here

  result = flex.double(flex.grid((Fmodel_indices.size(),500)))

  GF_FE1 = gen_fmodel_with_complex.from_structure(C2_structures[2],energy
           ).from_parameters(algorithm="direct")
  GF_FE1 = GF_FE1.as_P1_primitive()
  GF_FE2 = gen_fmodel_with_complex.from_structure(C2_structures[3],energy
           ).from_parameters(algorithm="direct")
  GF_FE2 = GF_FE2.as_P1_primitive()

  for incr in range(100):
    print ("incr is",incr)
    energy = 7070.5 + incr

    GF_FE1.reset_specific_at_energy(label_has="FE1",tables=FE1_model,newvalue=energy)
    # not sure if I need this, takes a lot of time # f_container = GF_FE1.get_fmodel()
    # not sure if I need this, takes a lot of time # Fcalc_FE1 = f_container.fmodel.f_calc()

    GF_FE2.reset_specific_at_energy(label_has="FE2",tables=FE2_model,newvalue=energy)
    # not sure if I need this, takes a lot of time # f_container = GF_FE2.get_fmodel()
    # not sure if I need this, takes a lot of time # Fcalc_FE2 = f_container.fmodel.f_calc()

    ALGO = structure_factors.from_scatterers(crystal_symmetry=CS,
                                             d_min=GF_whole7070.params2.high_resolution)
    #hack cctbx/xray/structure_factors/structure_factors_direct.h
    """
/*#if !defined(CCTBX_XRAY_STRUCTURE_FACTORS_DIRECT_NO_PRAGMA_OMP)
#if !defined(__DECCXX_VER) || (defined(_OPENMP) && _OPENMP > 199819)
        #pragma omp parallel for schedule(static)
#endif
#endif
*/
        #pragma omp parallel for
    """
    from_scatterers_direct_fe1 = ALGO(xray_structure=GF_FE1.xray_structure,
                                      miller_set=MS,algorithm="direct")
    Fcalc_FE1_dir = from_scatterers_direct_fe1.f_calc().data()
    from_scatterers_direct_fe2 = ALGO(xray_structure=GF_FE2.xray_structure,
                                      miller_set=MS,algorithm="direct")
    Fcalc_FE2_dir = from_scatterers_direct_fe2.f_calc().data()

    # Get total Fcalc at energy:
    F_bulk_non_Fe = base.matrix_copy_block(i_row=0,i_column=incr,
                    n_rows=Fmodel_indices.size(),n_columns=1)
    Fcalc_FE1_dir.reshape(flex.grid((Fmodel_indices.size(),1)))
    Fcalc_FE2_dir.reshape(flex.grid((Fmodel_indices.size(),1)))
    Fcalc_total = F_bulk_non_Fe + Fcalc_FE1_dir + Fcalc_FE2_dir

    result.matrix_paste_block_in_place(flex.norm(Fcalc_total),0,incr) # gives I = F * F

    gradient_flags=xray.structure_factors.gradient_flags(
      site=False,
      u_iso=False,
      u_aniso=False,
      occupancy=False,
      fp=True,
      fdp=True)
    xray.set_scatterer_grad_flags(scatterers = GF_FE1.xray_structure.scatterers(),
                                site       = gradient_flags.site,
                                u_iso      = gradient_flags.u_iso,
                                u_aniso    = gradient_flags.u_aniso,
                                occupancy  = gradient_flags.occupancy,
                                fp         = gradient_flags.fp,
                                fdp        = gradient_flags.fdp)
    xray.set_scatterer_grad_flags(scatterers = GF_FE2.xray_structure.scatterers(),
                                site       = gradient_flags.site,
                                u_iso      = gradient_flags.u_iso,
                                u_aniso    = gradient_flags.u_aniso,
                                occupancy  = gradient_flags.occupancy,
                                fp         = gradient_flags.fp,
                                fdp        = gradient_flags.fdp)

    #hack cctbx/xray/structure_factors/each_hkl_gradients_direct.h
    #pragma omp parallel for (line 226)
    sf1 = xray.ext.each_hkl_gradients_direct(
    MS.unit_cell(), MS.space_group(), MS.indices(), GF_FE1.xray_structure.scatterers(), None,
    GF_FE1.xray_structure.scattering_type_registry(), GF_FE1.xray_structure.site_symmetry_table(),
    0)
    sf2 = xray.ext.each_hkl_gradients_direct(
    MS.unit_cell(), MS.space_group(), MS.indices(), GF_FE2.xray_structure.scatterers(), None,
    GF_FE2.xray_structure.scattering_type_registry(), GF_FE2.xray_structure.site_symmetry_table(),
    0)

    sf1.d_fcalc_d_fp().reshape(flex.grid((Fmodel_indices.size(),1)))
    sf1.d_fcalc_d_fdp().reshape(flex.grid((Fmodel_indices.size(),1)))
    sf2.d_fcalc_d_fp().reshape(flex.grid((Fmodel_indices.size(),1)))
    sf2.d_fcalc_d_fdp().reshape(flex.grid((Fmodel_indices.size(),1)))
    parts = Fcalc_total.parts()
    A = parts[0]
    B = parts[1]
    for offset,vec2 in zip([100,200,300,400],
      [sf1.d_fcalc_d_fp(),sf1.d_fcalc_d_fdp(),sf2.d_fcalc_d_fp(),sf2.d_fcalc_d_fdp()]):
      vparts = vec2.parts()
      vA = vparts[0]
      vB = vparts[1]
      partial_I_partial_q = 2. * (A * vA + B * vB)
      result.matrix_paste_block_in_place(partial_I_partial_q,0,incr + offset)

  # result holds a table of structure factor intensities.  Rows are Miller indices H.
  # First 100 columns are I_H(energy, 100 channels). The next four groups of 100 columns
  # are the partial derivatives of I_H with respect to fp(FE1), fdp(FE1), fp(FE2), and fdp(FE2)
  # respectively all as a function of energy channel. Thus this is
  # the energy-dependent portion of the calculation that is dependent on the iron model.
  return result