atomes 1.1.16
atomes: an atomic scale modeling tool box
Loading...
Searching...
No Matches
rings-guttman.F90
Go to the documentation of this file.
1! This file is part of the 'atomes' software.
2!
3! 'atomes' is free software: you can redistribute it and/or modify it under the terms
4! of the GNU Affero General Public License as published by the Free Software Foundation,
5! either version 3 of the License, or (at your option) any later version.
6!
7! 'atomes' is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
8! without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
9! See the GNU General Public License for more details.
10!
11! You should have received a copy of the GNU Affero General Public License along with 'atomes'.
12! If not, see <https://www.gnu.org/licenses/>
13!
14! Copyright (C) 2022-2024 by CNRS and University of Strasbourg
15!
20
21INTEGER FUNCTION guttman_rings()
22
23USE parameters
24
25#ifdef OPENMP
26!$ USE OMP_LIB
27#endif
28IMPLICIT NONE
29
30#ifdef OPENMP
31INTEGER :: numth
32LOGICAL :: doatoms
33#endif
34
35INTERFACE
36 INTEGER FUNCTION recrings(VID)
37 INTEGER, INTENT(IN) :: vid
38 END FUNCTION
39END INTERFACE
40
42
43#ifdef OPENMP
44numth = omp_get_max_threads()
45doatoms=.false.
46! Dynamic allocation of pointers and tables used in the "RINGS" subroutine
47if (ns.ge.1 .and. ns.lt.numth) then
48 if (numth .ge. 2*(ns-1)) then
49 doatoms=.true.
50 else
51 numth=ns
52 endif
53endif
54
55if (all_atoms) doatoms=.true.
56
57if (doatoms) then
58 if (na.lt.numth) numth=na
59#ifdef DEBUG
60 write (6, *) "OpenMP on atoms, NUMTH= ",numth
61#endif
62 call guttman_ring_search_atoms (numth)
63else
64#ifdef DEBUG
65 write (6, *) "OpenMP on MD steps, NUMTH= ",numth
66#endif
67 call guttman_ring_search_steps (numth)
68endif
69#else
71#endif
72
74
75END FUNCTION
76
77#ifdef OPENMP
78SUBROUTINE guttman_ring_search_atoms (NUMTH)
79
80USE parameters
81
82!$ USE OMP_LIB
83IMPLICIT NONE
84
85INTEGER, INTENT(IN) :: numth
86TYPE (ring), DIMENSION(:), ALLOCATABLE :: the_ring
87INTEGER, DIMENSION(:), ALLOCATABLE :: tring, indte, indth
88INTEGER, DIMENSION(:,:,:), ALLOCATABLE :: savring, ordring
89INTEGER, DIMENSION(:,:,:), ALLOCATABLE :: savr, ordr
90INTEGER :: lora, lorb, lorc, ri
91
92INTERFACE
93 RECURSIVE SUBROUTINE inside_ring (THE_RING, FND, S_IR, AI_IR, RID, TAE, TAH, LRA, LRB, &
94 NRPAT, RSAVED, OSAVED, TRING, INDE, INDH, RESL, CPT, VPT)
95 USE parameters
96 TYPE (ring), DIMENSION(TAILLD), INTENT(INOUT) :: the_ring
97 LOGICAL, INTENT(INOUT) :: fnd
98 INTEGER, INTENT(IN) :: s_ir, ai_ir, rid
99 INTEGER, INTENT(INOUT) :: tae, tah
100 INTEGER, INTENT(IN) :: lra, lrb
101 INTEGER, DIMENSION(NA), INTENT(INOUT) :: nrpat
102 INTEGER, DIMENSION(TAILLR,NUMA,TAILLR), INTENT(INOUT) :: rsaved, osaved
103 INTEGER, DIMENSION(TAILLR), INTENT(INOUT) :: tring
104 INTEGER, DIMENSION(NUMA), INTENT(INOUT) :: inde, indh
105 INTEGER, DIMENSION(TAILLR), INTENT(INOUT) :: resl
106 INTEGER, DIMENSION(NA), INTENT(IN):: cpt
107 INTEGER, DIMENSION(NA,MAXN), INTENT(IN) :: vpt
108 END SUBROUTINE
109 INTEGER FUNCTION rings_to_ogl_menu (IDSEARCH, NRI)
110 USE parameters
111 INTEGER, INTENT(IN) :: idsearch
112 INTEGER, DIMENSION(TAILLR, NS), INTENT(IN) :: nri
113 END FUNCTION
114 INTEGER FUNCTION rings_to_ogl (STEP, IDSEARCH, NRI, RSAVED, OSAVED)
115 USE parameters
116 INTEGER, INTENT(IN) :: step, idsearch
117 INTEGER, DIMENSION(TAILLR,NS), INTENT(IN) :: nri
118 INTEGER, DIMENSION(TAILLR,NUMA,TAILLR), INTENT(IN) :: rsaved, osaved
119 END FUNCTION
120 SUBROUTINE del_this_ring (TLED, RSAVED, OSAVED, TRING, RESL, INDT)
121 USE parameters
122 INTEGER, INTENT(IN) :: tled
123 INTEGER, DIMENSION(TAILLR,NUMA,TAILLR), INTENT(INOUT) :: rsaved, osaved
124 INTEGER, DIMENSION(NUMA), INTENT(INOUT) :: indt
125 INTEGER, DIMENSION(TAILLR), INTENT(INOUT) :: resl, tring
126 END SUBROUTINE
127END INTERFACE
128
129ri = 0
130if(allocated(savring)) deallocate(savring)
131allocate(savring(taillr,numa,taillr), stat=err)
132if (err .ne. 0) then
133 alc_tab="SAVRING"
134 alc=.true.
135 goto 001
136endif
137if(allocated(ordring)) deallocate(ordring)
138allocate(ordring(taillr,numa,taillr), stat=err)
139if (err .ne. 0) then
140 alc_tab="ORDRING"
141 alc=.true.
142 goto 001
143endif
144if(allocated(cpat)) deallocate(cpat)
145allocate(cpat(na), stat=err)
146if (err .ne. 0) then
147 alc_tab="CPAT"
148 alc=.true.
149 goto 001
150endif
151 if(allocated(vpat)) deallocate(vpat)
152allocate(vpat(na,maxn), stat=err)
153if (err .ne. 0) then
154 alc_tab="VPAT"
155 alc=.true.
156 goto 001
157endif
158
159do i=1, ns
160
161 p=0
162
163 savring(:,:,:)=0
164 ordring(:,:,:)=0
165 call setup_cpat_vpat_ring (na, i, contj, voisj, cpat, vpat)
166
167 ! OpenMP on atoms only
168 !$OMP PARALLEL NUM_THREADS(NUMTH) DEFAULT (NONE) &
169 !$OMP& PRIVATE(MAXAT, MINAT, RPAT, APNA, SAUT, RUNSEARCH, &
170 !$OMP& j, k, l, m, n, o, LORA, LORB, LORC, TAILLE, TAILLH, THE_RING, RES_LIST, &
171 !$OMP& FOUND, ERR, SAVR, ORDR, TRING, INDTE, INDTH) &
172 !$OMP& SHARED(i, p, NUMTH, NS, NA, TLT, NSP, LOT, TAILLR, TAILLD, CONTJ, VOISJ, CPAT, VPAT, &
173 !$OMP& NUMA, FACTATRING, ATRING, MAXPNA, MINPNA, AMPAT, ABAB, NO_HOMO, ALLRINGS, &
174 !$OMP& TBR, ALC, ALC_TAB, NCELLS, PBC, MAXN, SAVRING, ORDRING, NRING, INDRING, PNA, ri)
175
176 if (allocated(rpat)) deallocate(rpat)
177 allocate(rpat(na), stat=err)
178 if (err .ne. 0) then
179 alc_tab="RPAT"
180 alc=.true.
181 goto 002
182 endif
183 if(allocated(res_list)) deallocate(res_list)
184 allocate(res_list(taillr), stat=err)
185 if (err .ne. 0) then
186 alc_tab="RES_LIST"
187 alc=.true.
188 goto 002
189 endif
190 if(allocated(indte)) deallocate(indte)
191 allocate(indte(numa), stat=err)
192 if (err .ne. 0) then
193 alc_tab="INDTE"
194 goto 002
195 endif
196 if(allocated(indth)) deallocate(indth)
197 allocate(indth(numa), stat=err)
198 if (err .ne. 0) then
199 alc_tab="INDTH"
200 goto 002
201 endif
202 if(allocated(apna)) deallocate(apna)
203 allocate(apna(taillr), stat=err)
204 if (err .ne. 0) then
205 alc_tab="APNA"
206 alc=.true.
207 goto 002
208 endif
209 if(allocated(savr)) deallocate(savr)
210 allocate(savr(taillr,numa,taillr), stat=err)
211 if (err .ne. 0) then
212 alc_tab="SAVR"
213 alc=.true.
214 goto 002
215 endif
216 if(allocated(ordr)) deallocate(ordr)
217 allocate(ordr(taillr,numa,taillr), stat=err)
218 if (err .ne. 0) then
219 alc_tab="ORDR"
220 alc=.true.
221 goto 002
222 endif
223 if(allocated(tring)) deallocate(tring)
224 allocate(tring(taillr), stat=err)
225 if (err .ne. 0) then
226 alc_tab="TRING"
227 alc=.true.
228 goto 002
229 endif
230 if(allocated(the_ring)) deallocate(the_ring)
231 allocate(the_ring(tailld), stat=err)
232 if (err .ne. 0) then
233 alc_tab="THE_RING"
234 alc=.true.
235 goto 002
236 endif
237
238 tring(:)=0
239 savr(:,:,:)=0
240 ordr(:,:,:)=0
241 !$OMP DO SCHEDULE(STATIC,NA/NUMTH)
242 do j=1, na
243
244 if (tbr .or. alc) goto 003
245 if (tlt .eq. nsp+1 .or. lot(j) .eq. tlt) then
246
247 !$OMP CRITICAL
248 p = p+1
249 o = p
250 !$OMP END CRITICAL
251 apna(:) = 0
252 minat=tailld
253 maxat=1
254 saut=.true.
255
256 if (cpat(j).ge.2) then
257
258 do l=1, cpat(j)
259
260 lora=lot(j)
261 lorb=lot(vpat(j,l))
262 if (abab) then
263 if (lora .ne. lorb) then
264 runsearch=.true.
265 else
266 runsearch=.false.
267 endif
268 else
269 runsearch=.true.
270 endif
271
272 if (runsearch) then
273
274 found=.false.
275 rpat(:) = 0
276 the_ring(1)%ATOM=j
277 the_ring(1)%SPEC=lora
278 the_ring(1)%NEIGHBOR=1
279 m = vpat(j,l)
280 rpat(m)=1
281 the_ring(2)%ATOM=m
282 the_ring(2)%SPEC=lorb
283 the_ring(2)%NEIGHBOR=cpat(m)
284 taille=taillr
285 taillh=taillr
286 res_list(:) = 0
287 indte(:) = 0
288 indth(:) = 0
289 call inside_ring (the_ring, found, i, m, 2, taille, taillh, lorb, lora, &
290 rpat, savr, ordr, tring, indte, indth, res_list, cpat, vpat)
291 if (alc) alc_tab="INSIDE_RING"
292 if (tbr .or. alc) goto 003
293 if (found) then
294 if (apna(taille) .eq. 0) then
295 apna(taille)=1
296 minat=min(minat,taille)
297 maxat=max(maxat,taille)
298 ! if (FACTATPNA) ATPNA(TAILLE,o,i)=1
299 saut=.false.
300 endif
301 elseif (no_homo .and. res_list(taillh) .ne. 0) then
302 ! Some shortest ring with HP bonds and < TAILLR were found
303 ! delete these rings now
304 call del_this_ring (taillh, savr, ordr, tring, res_list, indth)
305 else
306 if (contj(j,i) .ge. 2 .and. contj(m,i) .ge. 2) then
307 !$OMP ATOMIC
308 ampat(o,i)=ampat(o,i)+1
309 endif
310 endif
311 endif
312 enddo
313
314 endif
315
316 if (.not. saut) then
317 do k=3, taillr
318 do l=3, taillr
319 if (apna(k).eq.1 .and. apna(l).eq.1) then
320 !$OMP ATOMIC
321 pna(k,l,i)=pna(k,l,i)+1
322 endif
323 enddo
324 enddo
325 !$OMP ATOMIC
326 maxpna(maxat,i)=maxpna(maxat,i)+1
327 !$OMP ATOMIC
328 minpna(minat,i)=minpna(minat,i)+1
329 endif
330 endif
331
332 !$OMP CRITICAL
333 do k=3, taillr
334 if (tring(k).gt.0) then
335 if (nring(k,i).gt.0) then
336 o = 0
337 do l=1, tring(k)
338 do m=1, nring(k,i)
339 saut=.true.
340 do n=1, k
341 if (savring(k,m,n) .ne. savr(k,l,n)) then
342 saut=.false.
343 exit
344 endif
345 enddo
346 if (saut) exit
347 enddo
348 if (.not.saut) then
349 o = o + 1
350 if (nring(k,i)+o .gt. numa) then
351 tbr=.true.
352 goto 004
353 endif
354 do m=1, k
355 savring(k,nring(k,i)+o,m) = savr(k,l,m)
356 ordring(k,nring(k,i)+o,m) = ordr(k,l,m)
357 enddo
358 endif
359 enddo
360 nring(k,i)=nring(k,i)+o
361 else
362 do l=1, tring(k)
363 do m=1, k
364 savring(k,l,m) = savr(k,l,m)
365 ordring(k,l,m) = ordr(k,l,m)
366 enddo
367 enddo
368 nring(k,i) = tring(k)
369 endif
370 endif
371 enddo
372
373 004 continue
374 !$OMP END CRITICAL
375
376 003 continue
377 enddo
378 !$OMP END DO NOWAIT
379
380 002 continue
381 if (allocated(rpat)) deallocate (rpat)
382 if (allocated(res_list)) deallocate (res_list)
383 if (allocated(indte)) deallocate (indte)
384 if (allocated(indth)) deallocate (indth)
385 if (allocated(apna)) deallocate (apna)
386 if (allocated(tring)) deallocate (tring)
387 if (allocated(savr)) deallocate (savr)
388 if (allocated(ordr)) deallocate (ordr)
389 if (allocated(the_ring)) deallocate (the_ring)
390 !$OMP END PARALLEL
391 if (alc .or. tbr) goto 001
392
393 !do j=3, TAILLR
394 ! write (6, '("s= ",i4,", j= ",i2,", nr(",i2,",",i4,")= ",i2)') i,j,j,i, NRING(j,i)
395 ! if (NRING(j,i) .gt. 0) then
396 ! do k=1, NRING(j,i)
397 ! write (6, *) " k= ",k,", R(o)= ",ORDRING(j,k,1:j)
398 ! enddo
399 ! endif
400 !enddo
401 ri = ri + rings_to_ogl(i, 2, nring, savring, ordring)
402
403enddo
404
405001 continue
406
407if (alc) then
408 call show_error ("Impossible to allocate memory"//char(0), &
409 "Subroutine: GUTTMAN_RING_SEARCH_ATOMS"//char(0), "Table: "//alc_tab(1:len_trim(alc_tab))//char(0))
410endif
411
412if (allocated(cpat)) deallocate (cpat)
413if (allocated(vpat)) deallocate (vpat)
414if (allocated(savring)) deallocate (savring)
415if (allocated(ordring)) deallocate (ordring)
416
417if (ri .eq. ns) ri = rings_to_ogl_menu(2, nring)
418
419END SUBROUTINE
420#endif
421
422#ifdef OPENMP
423SUBROUTINE guttman_ring_search_steps (NUMTH)
424
425USE parameters
426!$ USE OMP_LIB
427IMPLICIT NONE
428INTEGER, INTENT(IN) :: numth
429#else
431
432USE parameters
433IMPLICIT NONE
434#endif
435TYPE (RING), DIMENSION(:), ALLOCATABLE :: THE_RING
436INTEGER, DIMENSION(:), ALLOCATABLE :: TRING, INDTE, INDTH
437INTEGER, DIMENSION(:,:,:), ALLOCATABLE :: SAVRING, ORDRING
438INTEGER :: LORA, LORB, ri
439
440INTERFACE
441 RECURSIVE SUBROUTINE inside_ring (THE_RING, FND, S_IR, AI_IR, RID, TAE, TAH, LRA, LRB, &
442 NRPAT, RSAVED, OSAVED, TRING, INDE, INDH, RESL, CPT, VPT)
443 USE parameters
444 TYPE (RING), DIMENSION(TAILLD), INTENT(INOUT) :: THE_RING
445 LOGICAL, INTENT(INOUT) :: FND
446 INTEGER, INTENT(IN) :: S_IR, AI_IR, RID
447 INTEGER, INTENT(INOUT) :: TAE, TAH
448 INTEGER, INTENT(IN) :: LRA, LRB
449 INTEGER, DIMENSION(NA), INTENT(INOUT) :: NRPAT
450 INTEGER, DIMENSION(TAILLR,NUMA,TAILLR), INTENT(INOUT) :: RSAVED, OSAVED
451 INTEGER, DIMENSION(TAILLR), INTENT(INOUT) :: TRING
452 INTEGER, DIMENSION(NUMA), INTENT(INOUT) :: INDE, INDH
453 INTEGER, DIMENSION(TAILLR), INTENT(INOUT) :: RESL
454 INTEGER, DIMENSION(NA), INTENT(IN):: CPT
455 INTEGER, DIMENSION(NA,MAXN), INTENT(IN) :: VPT
456 END SUBROUTINE
457 INTEGER FUNCTION rings_to_ogl_menu (IDSEARCH, NRI)
458 USE parameters
459 INTEGER, INTENT(IN) :: IDSEARCH
460 INTEGER, DIMENSION(TAILLR, NS), INTENT(IN) :: NRI
461 END FUNCTION
462 INTEGER FUNCTION rings_to_ogl (STEP, IDSEARCH, NRI, RSAVED, OSAVED)
463 USE parameters
464 INTEGER, INTENT(IN) :: STEP, IDSEARCH
465 INTEGER, DIMENSION(TAILLR,NS), INTENT(IN) :: NRI
466 INTEGER, DIMENSION(TAILLR,NUMA,TAILLR), INTENT(IN) :: RSAVED, OSAVED
467 END FUNCTION
468 SUBROUTINE del_this_ring (TLED, RSAVED, OSAVED, TRING, RESL, INDT)
469 USE parameters
470 INTEGER, INTENT(IN) :: TLED
471 INTEGER, DIMENSION(TAILLR,NUMA,TAILLR), INTENT(INOUT) :: RSAVED, OSAVED
472 INTEGER, DIMENSION(NUMA), INTENT(INOUT) :: INDT
473 INTEGER, DIMENSION(TAILLR), INTENT(INOUT) :: RESL, TRING
474 END SUBROUTINE
475END INTERFACE
476
477ri = 0
478
479#ifdef OPENMP
480! OpenMP on steps only
481!$OMP PARALLEL NUM_THREADS(NUMTH) DEFAULT (NONE) &
482!$OMP& PRIVATE(TAILLE, TAILLH, MAXAT, MINAT, SAUT, RUNSEARCH, &
483!$OMP& i, j, k, l, m, n, o, LORA, LORB, THE_RING, RES_LIST, INDTE, INDTH, APNA, &
484!$OMP& FOUND, ERR, TRING, SAVRING, ORDRING, RPAT, CPAT, VPAT) &
485!$OMP& SHARED(p, NUMTH, NS, NA, TLT, NSP, LOT, TAILLR, TAILLD, CONTJ, VOISJ, &
486!$OMP& NUMA, FACTATRING, ATRING, MAXPNA, MINPNA, AMPAT, ABAB, NO_HOMO, ALLRINGS, &
487!$OMP& TBR, ALC, ALC_TAB, NCELLS, THE_BOX, FULLPOS, PBC, MAXN, NRING, INDRING, PNA, ri)
488#endif
489
490if(allocated(savring)) deallocate(savring)
491allocate(savring(taillr,numa,taillr), stat=err)
492if (err .ne. 0) then
493 alc_tab="SAVRING"
494 alc=.true.
495 goto 001
496endif
497if(allocated(ordring)) deallocate(ordring)
498allocate(ordring(taillr,numa,taillr), stat=err)
499if (err .ne. 0) then
500 alc_tab="ORDRING"
501 alc=.true.
502 goto 001
503endif
504if(allocated(cpat)) deallocate(cpat)
505allocate(cpat(na), stat=err)
506if (err .ne. 0) then
507 alc_tab="CPAT"
508 alc=.true.
509 goto 001
510endif
511if(allocated(vpat)) deallocate(vpat)
512allocate(vpat(na,maxn), stat=err)
513if (err .ne. 0) then
514 alc_tab="VPAT"
515 alc=.true.
516 goto 001
517endif
518if (allocated(rpat)) deallocate(rpat)
519allocate(rpat(na), stat=err)
520if (err .ne. 0) then
521 alc_tab="RPAT"
522 alc=.true.
523 goto 001
524endif
525if(allocated(res_list)) deallocate(res_list)
526allocate(res_list(taillr), stat=err)
527if (err .ne. 0) then
528 alc_tab="RES_LIST"
529 alc=.true.
530 goto 001
531endif
532if(allocated(indte)) deallocate(indte)
533allocate(indte(numa), stat=err)
534if (err .ne. 0) then
535 alc_tab="INDTE"
536 goto 001
537endif
538if(allocated(indth)) deallocate(indth)
539allocate(indth(numa), stat=err)
540if (err .ne. 0) then
541 alc_tab="INDTH"
542 goto 001
543endif
544if(allocated(apna)) deallocate(apna)
545allocate(apna(taillr), stat=err)
546if (err .ne. 0) then
547 alc_tab="APNA"
548 alc=.true.
549 goto 001
550endif
551if(allocated(tring)) deallocate(tring)
552allocate(tring(taillr), stat=err)
553if (err .ne. 0) then
554 alc_tab="TRING"
555 alc=.true.
556 goto 001
557endif
558if(allocated(the_ring)) deallocate(the_ring)
559allocate(the_ring(tailld), stat=err)
560if (err .ne. 0) then
561 alc_tab="THE_RING"
562 alc=.true.
563 goto 001
564endif
565
566#ifdef OPENMP
567!$OMP DO SCHEDULE(STATIC,NS/NUMTH)
568#endif
569do i=1, ns
570
571 if (tbr .or. alc) goto 002
572 savring(:,:,:)=0
573 ordring(:,:,:)=0
574 tring(:)=0
575 call setup_cpat_vpat_ring (na, i, contj, voisj, cpat, vpat)
576
577 o=0
578 do j=1, na
579
580 if (tlt .eq. nsp+1 .or. lot(j) .eq. tlt) then
581
582 o=o+1
583 apna(:) = 0
584 minat=taillr
585 maxat=1
586 saut=.true.
587
588 if (cpat(j).ge.2) then
589
590 do l=1, cpat(j)
591
592 lora=lot(j)
593 lorb=lot(vpat(j,l))
594 if (abab) then
595 if (lora .ne. lorb) then
596 runsearch=.true.
597 else
598 runsearch=.false.
599 endif
600 else
601 runsearch=.true.
602 endif
603
604 if (runsearch) then
605 found=.false.
606 rpat(:) = 0
607 the_ring(1)%ATOM=j
608 the_ring(1)%SPEC=lora
609 the_ring(1)%NEIGHBOR=1
610 m = vpat(j,l)
611 rpat(m)=1
612 the_ring(2)%ATOM=m
613 the_ring(2)%SPEC=lorb
614 the_ring(2)%NEIGHBOR=cpat(m)
615 taille=taillr
616 taillh=taillr
617 res_list(:) = 0
618 indte(:) = 0
619 indth(:) = 0
620 call inside_ring (the_ring, found, i, m, 2, taille, taillh, lorb, lora, &
621 rpat, savring, ordring, tring, indte, indth, res_list, cpat, vpat)
622 if (alc) alc_tab="INSIDE_RING"
623 if (tbr .or. alc) goto 002
624 if (found) then
625 if (apna(taille) .eq. 0) then
626 apna(taille)=1
627 minat=min(minat,taille)
628 maxat=max(maxat,taille)
629 ! if (FACTATPNA) ATPNA(TAILLE,o,i)=1
630 saut=.false.
631 endif
632 elseif (no_homo .and. res_list(taillh) .ne. 0) then
633 ! Some shortest ring with HP bonds and < TAILLR were found
634 ! delete these rings now
635 call del_this_ring (taillh, savring, ordring, tring, res_list, indth)
636 else
637 if (contj(j,i) .ge. 2 .and. contj(m,i) .ge. 2) ampat(o,i)=ampat(o,i)+1
638 endif
639 endif
640 enddo
641
642 endif
643
644 if (.not. saut) then
645 do k=3, taillr
646 do l=3, taillr
647 if (apna(k).eq.1 .and. apna(l).eq.1) then
648 pna(k,l,i)=pna(k,l,i)+1
649 endif
650 enddo
651 enddo
652 maxpna(maxat,i)=maxpna(maxat,i)+1
653 minpna(minat,i)=minpna(minat,i)+1
654 endif
655 endif
656
657 enddo
658
659 do j=3, taillr
660 nring(j,i) = tring(j)
661 enddo
662
663 !do j=3, TAILLR
664 ! write (6, '("s= ",i4,", j= ",i2,", nr(",i2,",",i4,")= ",i2)') i,j,j,i, NRING(j,i)
665 ! if (NRING(j,i) .gt. 0) then
666 ! do k=1, NRING(j,i)
667 ! write (6, *) " k= ",k,", R(o)= ",ORDRING(j,k,1:j)
668 ! enddo
669 ! endif
670 !enddo
671 j = rings_to_ogl(i, 2, nring, savring, ordring)
672#ifdef OPENMP
673 !$OMP ATOMIC
674#endif
675 ri = ri + j
676 002 continue
677enddo
678#ifdef OPENMP
679!$OMP END DO NOWAIT
680#endif
681
682001 continue
683
684if (alc) then
685 call show_error ("Impossible to allocate memory"//char(0), &
686 "Subroutine: GUTTMAN_RING_SEARCH_STEPS"//char(0), "Table: "//alc_tab(1:len_trim(alc_tab))//char(0))
687endif
688
689if (allocated(tring)) deallocate (tring)
690if (allocated(the_ring)) deallocate (the_ring)
691if (allocated(savring)) deallocate (savring)
692if (allocated(ordring)) deallocate (ordring)
693if (allocated(cpat)) deallocate (cpat)
694if (allocated(vpat)) deallocate (vpat)
695if (allocated(rpat)) deallocate (rpat)
696if (allocated(res_list)) deallocate (res_list)
697if (allocated(indte)) deallocate (indte)
698if (allocated(indth)) deallocate (indth)
699if (allocated(apna)) deallocate (apna)
700
701#ifdef OPENMP
702!$OMP END PARALLEL
703#endif
704
705if (ri .eq. ns) ri = rings_to_ogl_menu(2, nring)
706
707END SUBROUTINE
#define min(a, b)
Definition global.h:81
#define max(a, b)
Definition global.h:80
void show_error(char *error, int val, GtkWidget *win)
show error message
Definition interface.c:293
integer function recrings(vid)
Definition resrings.F90:22
integer function guttman_rings()
subroutine guttman_ring_search_steps()
subroutine del_this_ring(tled, rsaved, osaved, tring, resl, indt)
subroutine setup_cpat_vpat_ring(nat, str, cont, vois, cpt, vpt)
recursive subroutine inside_ring(the_ring, fnd, s_ir, ai_ir, rid, tae, tah, lra, lrb, nrpat, rsaved, osaved, tring, inde, indh, resl, cpt, vpt)
integer function rings_to_ogl(step, idsearch, nri, rsaved, osaved)
Definition rings_ogl.F90:22
integer function rings_to_ogl_menu(idsearch, nri)