cluster & wham Adam's PBC corrections
[unres.git] / source / cluster / wham / src-HCD-5D / read_constr_homology.F
1       subroutine read_constr_homology
2
3       include 'DIMENSIONS'
4 #ifdef MPI
5       include 'mpif.h'
6 #endif
7       include 'COMMON.SETUP'
8       include 'COMMON.CONTROL'
9       include 'COMMON.CHAIN'
10       include 'COMMON.IOUNITS'
11       include 'COMMON.GEO'
12       include 'COMMON.INTERACT'
13       include 'COMMON.NAMES'
14       include 'COMMON.HOMRESTR'
15       include 'COMMON.HOMOLOGY'
16 c
17 c For new homol impl
18 c
19       include 'COMMON.VAR'
20 c     include 'include_unres/COMMON.VAR'
21 c
22
23 c     double precision odl_temp,sigma_odl_temp,waga_theta,waga_d,
24 c    &                 dist_cut
25 c     common /przechowalnia/ odl_temp(maxres,maxres,max_template),
26 c    &    sigma_odl_temp(maxres,maxres,max_template)
27       character*2 kic2
28       character*24 model_ki_dist, model_ki_angle
29       character*500 controlcard
30       integer ki, i, j, k, l, ii_in_use(maxdim),i_tmp,idomain_tmp
31       integer idomain(max_template,maxres)
32       logical lprn /.true./
33       integer ilen
34       external ilen
35       logical liiflag
36 c
37 c     FP - Nov. 2014 Temporary specifications for new vars
38 c
39       double precision rescore_tmp,x12,y12,z12,rescore2_tmp
40      &    rescore3_tmp
41       double precision, dimension (max_template,maxres) :: rescore
42       double precision, dimension (max_template,maxres) :: rescore2
43       double precision, dimension (max_template,maxres) :: rescore3
44       character*24 tpl_k_rescore
45 c -----------------------------------------------------------------
46 c Reading multiple PDB ref structures and calculation of retraints
47 c not using pre-computed ones stored in files model_ki_{dist,angle}
48 c FP (Nov., 2014)
49 c -----------------------------------------------------------------
50 c
51 c
52 c Alternative: reading from input
53       call card_concat(controlcard)
54       call reada(controlcard,"HOMOL_DIST",waga_dist,1.0d0)
55       call reada(controlcard,"HOMOL_ANGLE",waga_angle,1.0d0)
56       call reada(controlcard,"HOMOL_THETA",waga_theta,1.0d0) ! new
57       call reada(controlcard,"HOMOL_SCD",waga_d,1.0d0) ! new
58       call reada(controlcard,'DIST_CUT',dist_cut,5.0d0) ! for diff ways of calc sigma
59       call reada(controlcard,'DIST2_CUT',dist2_cut,9999.0d0)
60       call readi(controlcard,"HOMOL_NSET",homol_nset,1)
61       read2sigma=(index(controlcard,'READ2SIGMA').gt.0)
62       call readi(controlcard,"IHSET",ihset,1)       
63       write (iout,*) "homol_nset ",homol_nset
64       if (homol_nset.gt.1)then
65          call card_concat(controlcard)
66          read(controlcard,*) (waga_homology(i),i=1,homol_nset) 
67 c         if(me.eq.king .or. .not. out1file .and. fg_rank.eq.0) then
68 c          write(iout,*) "iset homology_weight "
69 c         do i=1,homol_nset
70 c          write(iout,*) i,waga_homology(i)
71 c         enddo
72 c         endif
73          iset=mod(kolor,homol_nset)+1
74       else
75        iset=1
76        waga_homology(1)=1.0
77       endif
78 c     write(iout,*) "waga_homology(",iset,")",waga_homology(iset)
79
80 cd      write (iout,*) "nnt",nnt," nct",nct
81 cd      call flush(iout)
82
83
84       lim_odl=0
85       lim_dih=0
86 c
87 c  New
88 c
89       lim_theta=0
90       lim_xx=0
91 c
92 c  Reading HM global scores (prob not required)
93 c
94       do i = nnt,nct
95         do k=1,constr_homology
96          idomain(k,i)=0
97         enddo
98       enddo
99 c     open (4,file="HMscore")
100 c     do k=1,constr_homology
101 c       read (4,*,end=521) hmscore_tmp
102 c       hmscore(k)=hmscore_tmp ! Another transformation can be used 
103 c       write(*,*) "Model", k, ":", hmscore(k)
104 c     enddo
105 c521  continue
106
107       ii=0
108       do i = nnt,nct-2 
109         do j=i+2,nct 
110         ii=ii+1
111         ii_in_use(ii)=0
112         enddo
113       enddo
114 c     write(iout,*) "waga_theta",waga_theta,"waga_d",waga_d
115
116       if (read_homol_frag) then
117         call read_klapaucjusz
118       else
119
120       do k=1,constr_homology
121
122         read(inp,'(a)') pdbfile
123 c  Next stament causes error upon compilation (?)
124 c       if(me.eq.king.or. .not. out1file)
125 c         write (iout,'(2a)') 'PDB data will be read from file ',
126 c    &   pdbfile(:ilen(pdbfile))
127          write (iout,'(a,5x,a)') 'HOMOL: Opening PDB file',
128      &  pdbfile(:ilen(pdbfile))
129         open(ipdbin,file=pdbfile,status='old',err=33)
130         goto 34
131   33    write (iout,'(a,5x,a)') 'Error opening PDB file',
132      &  pdbfile(:ilen(pdbfile))
133         stop
134   34    continue
135 c        print *,'Begin reading pdb data'
136 c
137 c Files containing res sim or local scores (former containing sigmas)
138 c
139
140         write(kic2,'(bz,i2.2)') k
141
142         tpl_k_rescore="template"//kic2//".sco"
143
144         unres_pdb=.false.
145         if (read2sigma) then
146           call readpdb_template(k)
147           close(ipdbin)
148         else
149           call readpdb(out_template_coord)
150           close(ipdbin)
151         endif
152
153 c        call readpdb
154         do i=1,2*nres
155           do j=1,3
156             crefjlee(j,i)=c(j,i)
157           enddo
158         enddo
159 #ifdef DEBUG
160         do i=1,nres
161           write (iout,'(i5,3f8.3,5x,3f8.3)') i,(crefjlee(j,i),j=1,3),
162      &      (crefjlee(j,i+nres),j=1,3)
163         enddo
164         write (iout,*) "read_constr_homology: after reading pdb file"
165         call flush(iout)
166 #endif
167
168 c
169 c     Distance restraints
170 c
171 c          ... --> odl(k,ii)
172 C Copy the coordinates from reference coordinates (?)
173         do i=1,2*nres
174           do j=1,3
175             c(j,i)=cref(j,i)
176 c           write (iout,*) "c(",j,i,") =",c(j,i)
177           enddo
178         enddo
179 c
180 c From read_dist_constr (commented out 25/11/2014 <-> res sim)
181 c
182 c         write(iout,*) "tpl_k_rescore - ",tpl_k_rescore
183           open (ientin,file=tpl_k_rescore,status='old')
184           if (nnt.gt.1) rescore(k,1)=0.0d0
185           do irec=nnt,nct ! loop for reading res sim 
186             if (read2sigma) then
187              read (ientin,*,end=1401) i_tmp,rescore2_tmp,rescore_tmp,
188      &                                rescore3_tmp,idomain_tmp
189              i_tmp=i_tmp+nnt-1
190              idomain(k,i_tmp)=idomain_tmp
191              rescore(k,i_tmp)=rescore_tmp
192              rescore2(k,i_tmp)=rescore2_tmp
193              rescore3(k,i_tmp)=rescore3_tmp
194              write(iout,'(a7,i5,3f10.5,i5)') "rescore",
195      &                      i_tmp,rescore2_tmp,rescore_tmp,
196      &                                rescore3_tmp,idomain_tmp
197             else
198              idomain(k,irec)=1
199              read (ientin,*,end=1401) rescore_tmp
200
201 c           rescore(k,irec)=rescore_tmp+1.0d0 ! to avoid 0 values
202              rescore(k,irec)=0.5d0*(rescore_tmp+0.5d0) ! alt transf to reduce scores
203 c           write(iout,*) "rescore(",k,irec,") =",rescore(k,irec)
204             endif
205           enddo  
206  1401   continue
207         close (ientin)        
208         if (waga_dist.ne.0.0d0) then
209           ii=0
210           do i = nnt,nct-2 
211             do j=i+2,nct 
212
213               x12=c(1,i)-c(1,j)
214               y12=c(2,i)-c(2,j)
215               z12=c(3,i)-c(3,j)
216               distal=dsqrt(x12*x12+y12*y12+z12*z12) 
217 c              write (iout,*) k,i,j,distal,dist2_cut
218
219             if (idomain(k,i).eq.idomain(k,j).and.idomain(k,i).ne.0
220      &            .and. distal.le.dist2_cut ) then
221
222               ii=ii+1
223               ii_in_use(ii)=1
224               l_homo(k,ii)=.true.
225
226 c             write (iout,*) "k",k
227 c             write (iout,*) "i",i," j",j," constr_homology",
228 c    &                       constr_homology
229               ires_homo(ii)=i
230               jres_homo(ii)=j
231               odl(k,ii)=distal
232               if (read2sigma) then
233                 sigma_odl(k,ii)=0
234                 do ik=i,j
235                  sigma_odl(k,ii)=sigma_odl(k,ii)+rescore2(k,ik)
236                 enddo
237                 sigma_odl(k,ii)=sigma_odl(k,ii)/(j-i+1)
238                 if (odl(k,ii).gt.dist_cut) sigma_odl(k,ii) = 
239      &        sigma_odl(k,ii)*dexp(0.5d0*(odl(k,ii)/dist_cut)**2-0.5d0)
240               else
241                 if (odl(k,ii).le.dist_cut) then
242                  sigma_odl(k,ii)=rescore(k,i)+rescore(k,j) 
243                 else
244 #ifdef OLDSIGMA
245                  sigma_odl(k,ii)=(rescore(k,i)+rescore(k,j))* 
246      &                      dexp(0.5d0*(odl(k,ii)/dist_cut)**2)
247 #else
248                  sigma_odl(k,ii)=(rescore(k,i)+rescore(k,j))* 
249      &                      dexp(0.5d0*(odl(k,ii)/dist_cut)**2-0.5d0)
250 #endif
251                 endif
252               endif
253               sigma_odl(k,ii)=1.0d0/(sigma_odl(k,ii)*sigma_odl(k,ii)) 
254             else
255               ii=ii+1
256               l_homo(k,ii)=.false.
257             endif
258             enddo
259           enddo
260         lim_odl=ii
261         endif
262 c
263 c     Theta, dihedral and SC retraints
264 c
265         if (waga_angle.gt.0.0d0) then
266 c         open (ientin,file=tpl_k_sigma_dih,status='old')
267 c         do irec=1,maxres-3 ! loop for reading sigma_dih
268 c            read (ientin,*,end=1402) i,j,ki,l,sigma_dih(k,i+nnt-1) ! j,ki,l what for?
269 c            if (i+nnt-1.gt.lim_dih) lim_dih=i+nnt-1 ! right?
270 c            sigma_dih(k,i+nnt-1)=sigma_dih(k,i+nnt-1)* ! not inverse because of use of res. similarity
271 c    &                            sigma_dih(k,i+nnt-1)
272 c         enddo
273 c1402   continue
274 c         close (ientin)
275           do i = nnt+3,nct 
276             if (idomain(k,i).eq.0) then 
277                sigma_dih(k,i)=0.0
278                cycle
279             endif
280             dih(k,i)=phiref(i) ! right?
281 c           read (ientin,*) sigma_dih(k,i) ! original variant
282 c             write (iout,*) "dih(",k,i,") =",dih(k,i)
283 c             write(iout,*) "rescore(",k,i,") =",rescore(k,i),
284 c    &                      "rescore(",k,i-1,") =",rescore(k,i-1),
285 c    &                      "rescore(",k,i-2,") =",rescore(k,i-2),
286 c    &                      "rescore(",k,i-3,") =",rescore(k,i-3)
287
288             sigma_dih(k,i)=(rescore(k,i)+rescore(k,i-1)+
289      &                     rescore(k,i-2)+rescore(k,i-3))/4.0
290 c            if (read2sigma) sigma_dih(k,i)=sigma_dih(k,i)/4.0
291 c           write (iout,*) "Raw sigmas for dihedral angle restraints"
292 c           write (iout,'(i5,10(2f8.2,4x))') i,sigma_dih(k,i)
293 c           sigma_dih(k,i)=hmscore(k)*rescore(k,i)*rescore(k,i-1)*
294 c                          rescore(k,i-2)*rescore(k,i-3)  !  right expression ?
295 c   Instead of res sim other local measure of b/b str reliability possible
296             if (sigma_dih(k,i).ne.0)
297      &      sigma_dih(k,i)=1.0d0/(sigma_dih(k,i)*sigma_dih(k,i))
298 c           sigma_dih(k,i)=sigma_dih(k,i)*sigma_dih(k,i)
299           enddo
300           lim_dih=nct-nnt-2 
301         endif
302
303         if (waga_theta.gt.0.0d0) then
304 c         open (ientin,file=tpl_k_sigma_theta,status='old')
305 c         do irec=1,maxres-2 ! loop for reading sigma_theta, right bounds?
306 c            read (ientin,*,end=1403) i,j,ki,sigma_theta(k,i+nnt-1) ! j,ki what for?
307 c            sigma_theta(k,i+nnt-1)=sigma_theta(k,i+nnt-1)* ! not inverse because of use of res. similarity
308 c    &                              sigma_theta(k,i+nnt-1)
309 c         enddo
310 c1403   continue
311 c         close (ientin)
312
313           do i = nnt+2,nct ! right? without parallel.
314 c         do i = i=1,nres ! alternative for bounds acc to readpdb?
315 c         do i=ithet_start,ithet_end ! with FG parallel.
316              if (idomain(k,i).eq.0) then  
317               sigma_theta(k,i)=0.0
318               cycle
319              endif
320              thetatpl(k,i)=thetaref(i)
321 c            write (iout,*) "thetatpl(",k,i,") =",thetatpl(k,i)
322 c            write(iout,*)  "rescore(",k,i,") =",rescore(k,i),
323 c    &                      "rescore(",k,i-1,") =",rescore(k,i-1),
324 c    &                      "rescore(",k,i-2,") =",rescore(k,i-2)
325 c            read (ientin,*) sigma_theta(k,i) ! 1st variant
326              sigma_theta(k,i)=(rescore(k,i)+rescore(k,i-1)+
327      &                        rescore(k,i-2))/3.0
328 c             if (read2sigma) sigma_theta(k,i)=sigma_theta(k,i)/3.0
329              if (sigma_theta(k,i).ne.0)
330      &       sigma_theta(k,i)=1.0d0/(sigma_theta(k,i)*sigma_theta(k,i))
331
332 c            sigma_theta(k,i)=hmscore(k)*rescore(k,i)*rescore(k,i-1)*
333 c                             rescore(k,i-2) !  right expression ?
334 c            sigma_theta(k,i)=sigma_theta(k,i)*sigma_theta(k,i)
335           enddo
336         endif
337
338         if (waga_d.gt.0.0d0) then
339 c       open (ientin,file=tpl_k_sigma_d,status='old')
340 c         do irec=1,maxres-1 ! loop for reading sigma_theta, right bounds?
341 c            read (ientin,*,end=1404) i,j,sigma_d(k,i+nnt-1) ! j,ki what for?
342 c            sigma_d(k,i+nnt-1)=sigma_d(k,i+nnt-1)* ! not inverse because of use of res. similarity
343 c    &                          sigma_d(k,i+nnt-1)
344 c         enddo
345 c1404   continue
346
347           do i = nnt,nct ! right? without parallel.
348 c         do i=2,nres-1 ! alternative for bounds acc to readpdb?
349 c         do i=loc_start,loc_end ! with FG parallel.
350                if (itype(i).eq.10) cycle 
351                if (idomain(k,i).eq.0 ) then 
352                   sigma_d(k,i)=0.0
353                   cycle
354                endif
355                xxtpl(k,i)=xxref(i)
356                yytpl(k,i)=yyref(i)
357                zztpl(k,i)=zzref(i)
358 c              write (iout,*) "xxtpl(",k,i,") =",xxtpl(k,i)
359 c              write (iout,*) "yytpl(",k,i,") =",yytpl(k,i)
360 c              write (iout,*) "zztpl(",k,i,") =",zztpl(k,i)
361 c              write(iout,*)  "rescore(",k,i,") =",rescore(k,i)
362 c               sigma_d(k,i)=rescore(k,i) !  right expression ?
363                sigma_d(k,i)=rescore3(k,i) !  right expression ?
364                if (sigma_d(k,i).ne.0)
365      &          sigma_d(k,i)=1.0d0/(sigma_d(k,i)*sigma_d(k,i))
366
367 c              sigma_d(k,i)=hmscore(k)*rescore(k,i) !  right expression ?
368 c              sigma_d(k,i)=sigma_d(k,i)*sigma_d(k,i)
369 c              read (ientin,*) sigma_d(k,i) ! 1st variant
370           enddo
371         endif
372       enddo
373 c
374 c remove distance restraints not used in any model from the list
375 c shift data in all arrays
376 c
377       if (waga_dist.ne.0.0d0) then
378         ii=0
379         liiflag=.true.
380         do i=nnt,nct-2 
381          do j=i+2,nct 
382           ii=ii+1
383           if (ii_in_use(ii).eq.0.and.liiflag) then
384             liiflag=.false.
385             iistart=ii
386           endif
387           if (ii_in_use(ii).ne.0.and..not.liiflag.or.
388      &                   .not.liiflag.and.ii.eq.lim_odl) then
389              if (ii.eq.lim_odl) then
390               iishift=ii-iistart+1
391              else
392               iishift=ii-iistart
393              endif
394              liiflag=.true.
395              do ki=iistart,lim_odl-iishift
396               ires_homo(ki)=ires_homo(ki+iishift)
397               jres_homo(ki)=jres_homo(ki+iishift)
398               ii_in_use(ki)=ii_in_use(ki+iishift)
399               do k=1,constr_homology
400                odl(k,ki)=odl(k,ki+iishift)
401                sigma_odl(k,ki)=sigma_odl(k,ki+iishift)
402                l_homo(k,ki)=l_homo(k,ki+iishift)
403               enddo
404              enddo
405              ii=ii-iishift
406              lim_odl=lim_odl-iishift
407           endif
408          enddo
409         enddo
410       endif
411
412       endif ! .not. klapaucjusz     
413
414       if (constr_homology.gt.0) call homology_partition
415       if (constr_homology.gt.0) call init_int_table
416 cd      write (iout,*) "homology_partition: lim_theta= ",lim_theta,
417 cd     & "lim_xx=",lim_xx
418 c     write (iout,*) "ithet_start =",ithet_start,"ithet_end =",ithet_end
419 c     write (iout,*) "loc_start =",loc_start,"loc_end =",loc_end
420 c
421 c Print restraints
422 c
423       if (.not.lprn) return
424 cd      write(iout,*) "waga_theta",waga_theta,"waga_d",waga_d
425 c      if(me.eq.king .or. .not. out1file .and. fg_rank.eq.0) then
426        write (iout,*) "Distance restraints from templates"
427        do ii=1,lim_odl
428        write(iout,'(3i5,100(2f8.2,1x,l1,4x))') 
429      &  ii,ires_homo(ii),jres_homo(ii),
430      &  (odl(ki,ii),1.0d0/dsqrt(sigma_odl(ki,ii)),l_homo(ki,ii),
431      &  ki=1,constr_homology)
432        enddo
433        write (iout,*) "Dihedral angle restraints from templates"
434        do i=nnt+3,nct
435         write (iout,'(i5,a4,100(2f8.2,4x))') i,restyp(itype(i)),
436      &      (rad2deg*dih(ki,i),
437      &      rad2deg/dsqrt(sigma_dih(ki,i)),ki=1,constr_homology)
438        enddo
439        write (iout,*) "Virtual-bond angle restraints from templates"
440        do i=nnt+2,nct
441         write (iout,'(i5,a4,100(2f8.2,4x))') i,restyp(itype(i)),
442      &      (rad2deg*thetatpl(ki,i),
443      &      rad2deg/dsqrt(sigma_theta(ki,i)),ki=1,constr_homology)
444        enddo
445        write (iout,*) "SC restraints from templates"
446        do i=nnt,nct
447         write(iout,'(i5,100(4f8.2,4x))') i,
448      &  (xxtpl(ki,i),yytpl(ki,i),zztpl(ki,i),
449      &   1.0d0/dsqrt(sigma_d(ki,i)),ki=1,constr_homology)
450        enddo
451 c      endif
452 c -----------------------------------------------------------------
453       return
454       end
455 c----------------------------------------------------------------------
456       subroutine read_klapaucjusz
457
458       include 'DIMENSIONS'
459 #ifdef MPI
460       include 'mpif.h'
461 #endif
462       include 'COMMON.SETUP'
463       include 'COMMON.CONTROL'
464       include 'COMMON.HOMOLOGY'
465       include 'COMMON.CHAIN'
466       include 'COMMON.IOUNITS'
467       include 'COMMON.GEO'
468       include 'COMMON.INTERACT'
469       include 'COMMON.NAMES'
470       include 'COMMON.HOMRESTR'
471       character*256 fragfile
472       integer ninclust(maxclust),inclust(max_template,maxclust),
473      &  nresclust(maxclust),iresclust(maxres,maxclust)
474
475       character*2 kic2
476       character*24 model_ki_dist, model_ki_angle
477       character*500 controlcard
478       integer ki, i, j, k, l, ii_in_use(maxdim),i_tmp,idomain_tmp
479       integer idomain(max_template,maxres)
480       logical lprn /.true./
481       integer ilen
482       external ilen
483       logical liiflag
484 c
485 c
486       double precision rescore_tmp,x12,y12,z12,rescore2_tmp
487       double precision, dimension (max_template,maxres) :: rescore
488       double precision, dimension (max_template,maxres) :: rescore2
489       character*24 tpl_k_rescore
490
491 c
492 c For new homol impl
493 c
494       include 'COMMON.VAR'
495 c
496       call getenv("FRAGFILE",fragfile) 
497       write (iout,*) "read_klapaucjusz ",fragfile
498       open(ientin,file=fragfile,status="old",err=10)
499       read(ientin,*) constr_homology,nclust
500       l_homo = .false.
501       sigma_theta=0.0
502       sigma_d=0.0
503       sigma_dih=0.0
504 c Read pdb files 
505       do k=1,constr_homology 
506         read(ientin,'(a)') pdbfile
507         write (iout,'(a,5x,a)') 'KLAPAUCJUSZ: Opening PDB file',
508      &  pdbfile(:ilen(pdbfile))
509         open(ipdbin,file=pdbfile,status='old',err=33)
510         goto 34
511   33    write (iout,'(a,5x,a)') 'Error opening PDB file',
512      &  pdbfile(:ilen(pdbfile))
513         stop
514   34    continue
515         unres_pdb=.false.
516         call readpdb_template(k)
517 c        do i=1,2*nres
518 c          do j=1,3
519 c            chomo(j,i,k)=c(j,i)
520 c          enddo
521 c        enddo
522         do i=1,nres
523           rescore(k,i)=0.2d0
524           rescore2(k,i)=1.0d0
525         enddo
526       enddo
527 c Read clusters
528       do i=1,nclust
529         read(ientin,*) ninclust(i),nresclust(i)
530         read(ientin,*) (inclust(k,i),k=1,ninclust(i))
531         read(ientin,*) (iresclust(k,i),k=1,nresclust(i))
532       enddo
533 c
534 c Loop over clusters
535 c
536       do l=1,nclust
537         do ll = 1,ninclust(l)
538         
539         k = inclust(ll,l)
540         do i=1,nres
541           idomain(k,i)=0
542         enddo
543         do i=1,nresclust(l)
544           if (nnt.gt.1)  then
545             idomain(k,iresclust(i,l)+1) = 1
546           else
547             idomain(k,iresclust(i,l)) = 1
548           endif
549         enddo
550 c
551 c     Distance restraints
552 c
553 c          ... --> odl(k,ii)
554 C Copy the coordinates from reference coordinates (?)
555         do i=1,2*nres
556           do j=1,3
557             c(j,i)=chomo(j,i,k)
558 c           write (iout,*) "c(",j,i,") =",c(j,i)
559           enddo
560         enddo
561         call int_from_cart(.true.,.false.)
562         call sc_loc_geom(.false.)
563         do i=1,nres
564           thetaref(i)=theta(i)
565           phiref(i)=phi(i)
566         enddo
567         if (waga_dist.ne.0.0d0) then
568           ii=0
569           do i = nnt,nct-2 
570             do j=i+2,nct 
571
572               x12=c(1,i)-c(1,j)
573               y12=c(2,i)-c(2,j)
574               z12=c(3,i)-c(3,j)
575               distal=dsqrt(x12*x12+y12*y12+z12*z12) 
576 c              write (iout,*) k,i,j,distal,dist2_cut
577
578             if (idomain(k,i).eq.idomain(k,j).and.idomain(k,i).ne.0
579      &            .and. distal.le.dist2_cut ) then
580
581               ii=ii+1
582               ii_in_use(ii)=1
583               l_homo(k,ii)=.true.
584
585 c             write (iout,*) "k",k
586 c             write (iout,*) "i",i," j",j," constr_homology",
587 c    &                       constr_homology
588               ires_homo(ii)=i
589               jres_homo(ii)=j
590               odl(k,ii)=distal
591               if (read2sigma) then
592                 sigma_odl(k,ii)=0
593                 do ik=i,j
594                  sigma_odl(k,ii)=sigma_odl(k,ii)+rescore2(k,ik)
595                 enddo
596                 sigma_odl(k,ii)=sigma_odl(k,ii)/(j-i+1)
597                 if (odl(k,ii).gt.dist_cut) sigma_odl(k,ii) = 
598      &        sigma_odl(k,ii)*dexp(0.5d0*(odl(k,ii)/dist_cut)**2-0.5d0)
599               else
600                 if (odl(k,ii).le.dist_cut) then
601                  sigma_odl(k,ii)=rescore(k,i)+rescore(k,j) 
602                 else
603 #ifdef OLDSIGMA
604                  sigma_odl(k,ii)=(rescore(k,i)+rescore(k,j))* 
605      &                      dexp(0.5d0*(odl(k,ii)/dist_cut)**2)
606 #else
607                  sigma_odl(k,ii)=(rescore(k,i)+rescore(k,j))* 
608      &                      dexp(0.5d0*(odl(k,ii)/dist_cut)**2-0.5d0)
609 #endif
610                 endif
611               endif
612               sigma_odl(k,ii)=1.0d0/(sigma_odl(k,ii)*sigma_odl(k,ii)) 
613             else
614               ii=ii+1
615 c              l_homo(k,ii)=.false.
616             endif
617             enddo
618           enddo
619         lim_odl=ii
620         endif
621 c
622 c     Theta, dihedral and SC retraints
623 c
624         if (waga_angle.gt.0.0d0) then
625           do i = nnt+3,nct 
626             if (idomain(k,i).eq.0) then 
627 c               sigma_dih(k,i)=0.0
628                cycle
629             endif
630             dih(k,i)=phiref(i)
631             sigma_dih(k,i)=(rescore(k,i)+rescore(k,i-1)+
632      &                     rescore(k,i-2)+rescore(k,i-3))/4.0
633 c            write (iout,*) "k",k," l",l," i",i," rescore",rescore(k,i),
634 c     &       " sigma_dihed",sigma_dih(k,i)
635             if (sigma_dih(k,i).ne.0)
636      &       sigma_dih(k,i)=1.0d0/(sigma_dih(k,i)*sigma_dih(k,i))
637           enddo
638           lim_dih=nct-nnt-2 
639         endif
640
641         if (waga_theta.gt.0.0d0) then
642           do i = nnt+2,nct
643              if (idomain(k,i).eq.0) then  
644 c              sigma_theta(k,i)=0.0
645               cycle
646              endif
647              thetatpl(k,i)=thetaref(i)
648              sigma_theta(k,i)=(rescore(k,i)+rescore(k,i-1)+
649      &                        rescore(k,i-2))/3.0
650              if (sigma_theta(k,i).ne.0)
651      &       sigma_theta(k,i)=1.0d0/(sigma_theta(k,i)*sigma_theta(k,i))
652           enddo
653         endif
654
655         if (waga_d.gt.0.0d0) then
656           do i = nnt,nct
657                if (itype(i).eq.10) cycle 
658                if (idomain(k,i).eq.0 ) then 
659 c                  sigma_d(k,i)=0.0
660                   cycle
661                endif
662                xxtpl(k,i)=xxref(i)
663                yytpl(k,i)=yyref(i)
664                zztpl(k,i)=zzref(i)
665                sigma_d(k,i)=rescore(k,i)
666                if (sigma_d(k,i).ne.0)
667      &          sigma_d(k,i)=1.0d0/(sigma_d(k,i)*sigma_d(k,i))
668                if (i-nnt+1.gt.lim_xx) lim_xx=i-nnt+1
669           enddo
670         endif
671       enddo ! l
672       enddo ! ll
673 c
674 c remove distance restraints not used in any model from the list
675 c shift data in all arrays
676 c
677       if (waga_dist.ne.0.0d0) then
678         ii=0
679         liiflag=.true.
680         do i=nnt,nct-2 
681          do j=i+2,nct 
682           ii=ii+1
683           if (ii_in_use(ii).eq.0.and.liiflag) then
684             liiflag=.false.
685             iistart=ii
686           endif
687           if (ii_in_use(ii).ne.0.and..not.liiflag.or.
688      &                   .not.liiflag.and.ii.eq.lim_odl) then
689              if (ii.eq.lim_odl) then
690               iishift=ii-iistart+1
691              else
692               iishift=ii-iistart
693              endif
694              liiflag=.true.
695              do ki=iistart,lim_odl-iishift
696               ires_homo(ki)=ires_homo(ki+iishift)
697               jres_homo(ki)=jres_homo(ki+iishift)
698               ii_in_use(ki)=ii_in_use(ki+iishift)
699               do k=1,constr_homology
700                odl(k,ki)=odl(k,ki+iishift)
701                sigma_odl(k,ki)=sigma_odl(k,ki+iishift)
702                l_homo(k,ki)=l_homo(k,ki+iishift)
703               enddo
704              enddo
705              ii=ii-iishift
706              lim_odl=lim_odl-iishift
707           endif
708          enddo
709         enddo
710       endif
711 #ifdef DEBUG
712       write (iout,*) "ires_homo and jres_homo arrays, lim_odl",lim_odl
713       do i=1,lim_odl
714         write (iout,*) i,ires_homo(i),jres_homo(i)
715       enddo
716 #endif
717       return
718    10 stop "Error in fragment file"
719       end