The TELEMAC-MASCARET system  trunk
pares3d.F
Go to the documentation of this file.
1 ! ******************
2  SUBROUTINE pares3d
3 ! ******************
4 
5  &(nameinp,namelog,nparts,pmethod,format_med)
6 !
7 !**********************************************************************
8 ! 12/11/2009 CHRISTOPHE DENIS SINETICS/I23
9 ! NEW VERSION TO DECREASE THE PARES3D COMPUTING TIME BY IMPROVING
10 !
11 ! - THE TETRA-TRIA CONNECTION
12 ! - THE POSTPROCESSING
13 !
14 ! COMMENTS ON THIS NEW VERSION -> CD
15 ! *********************************************************************
16 !***********************************************************************
17 ! PARTEL VERSION 5.6 08/06/06 O.BOITEAU/F.DECUNG(SINETICS/LNHE)
18 ! PARTEL VERSION 5.8 02/07/07 F.DECUNG(LNHE)
19 ! VERSION DE DEVELOPPEMENT POUR PRISE EN COMPTE PB DECOUPAGE
20 ! F.DECUNG/O.BOITEAU (JANV 2008)
21 ! AJOUT LECTURE DU FORMAT MED (SEPT 2013) V.STOBIAC
22 ! COPYRIGHT 2006
23 !***********************************************************************
24 !
25 ! CONSTRUCTIONS DES FICHIERS POUR ALIMENTER LE FLOT DE DONNEES
26 ! PARALLELE LORS D'UN CALCUL ESTEL3D PARALLELE EN ECOULEMENT
27 !
28 !
29 !-----------------------------------------------------------------------
30 ! ARGUMENTS
31 ! .________________.____.______________________________________________.
32 ! | NOM |MODE| ROLE |
33 ! |________________|____|______________________________________________|
34 ! | NAMEINP | -->| NOM DU FICHIER DE GEOMETRIE ESTEL3D
35 ! | NAMELOG | -->| NOM DU FICHIER LOG
36 ! | NPARTS | -->| NOMBRE DE PARTITION
37 ! | PMETHOD | -->| METHODE DE PARTITIONEMENT 1 metis 2 scotch
38 ! | FORMAT_MED | -->| BOULEEN POUR LE MAILLAGE AU FORMAT MED
39 ! |________________|____|______________________________________________|
40 ! MODE: -->(DONNEE NON MODIFIEE),<--(RESULTAT),<-->(DONNEE MODIFIEE)
41 !
42 !-----------------------------------------------------------------------
43 !
44 ! APPELE PAR : PARTEL
45 !
46 ! SOUS-PROGRAMME APPELE :
47 ! PARTEL_ALLOER, PARTEL_ALLOER2 (GESTION MSGS)
48 ! METIS_PARTMESHDUAL (FROM METIS LIBRARY)
49 !***********************************************************************
50 !
52  USE util_pares
53 !
55  IMPLICIT NONE
56 !
57 #if defined HAVE_MED
58  include 'med.hf'
59 #endif
60 !
61 !
62 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
63 !
64  CHARACTER(LEN=MAXLENHARD), INTENT(IN) :: NAMEINP
65  CHARACTER(LEN=MAXLENHARD), INTENT(IN) :: NAMELOG
66  INTEGER, INTENT(IN) :: NPARTS
67  INTEGER, INTENT(IN) :: PMETHOD
68  LOGICAL, INTENT(IN) :: FORMAT_MED
69 !
70 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
71 !
72 !
73 ! VARIABLES LOCALES
74  CHARACTER(LEN=MAXLENHARD) :: NAMEINP2,NAMELOG2
75  INTEGER :: NINP=10,nlog=11,ninp2=12,nlog2=13
76  INTEGER :: I,I_LENINP,IERR,J,K,COMPT,
77  & n,numtet,numtri,numtrig,i_lenlog,l,ni,nf,nt,ibid,idd,
78  & compt1,compt2,compt3,nbtriidd,m,color1,
79  & color2,pr1,pr2,nbtetj,iddnt,nit,nft,mt,
80  & numtrib,numtetb,ibidc,nbretouche,indpu(1)
81  LOGICAL :: IS,LINTER
82  CHARACTER(LEN=300) :: TEXTERROR ! TEXTE MSG D'ERREUR
83  CHARACTER(LEN=8) :: STR8 ! TEXTE MSG D'ERREUR
84  CHARACTER(LEN=300) :: STR26 ! TEXTE MSG D'ERREUR
85  CHARACTER(LEN=80) :: TITRE ! MESH TITLE IN THE FILE
86  CHARACTER(LEN=2) :: MOINS1 ! "-1"
87  CHARACTER(LEN=4) :: BLANC ! WHITE SPACE
88 
89  ! ADDITION JP RENAUD 15/02/2007
90  CHARACTER(LEN=200) :: LINE ! ONE LINE, 200 CHARACTERS MAXADDCH
91  INTEGER :: POS ! POSITION OF A CHARACTER IN THE LINE
92  INTEGER :: IOS ! STATUS INTEGER
93  ! END ADDITION JP RENAUD
94  CHARACTER(LEN=72) :: THEFORMAT
95 
96  CHARACTER(LEN=80), ALLOCATABLE :: LOGFAMILY(:) ! LOG INFORMATIONS
97  INTEGER :: NSEC ! TYPE OF THE SECTION READ
98  INTEGER, PARAMETER :: NSEC1=151 ! MESH TITLE SECTION ID
99  INTEGER, PARAMETER :: NSEC2=2411 ! NODES COORDINATES SECTION ID
100  INTEGER, PARAMETER :: NSEC3=2412 ! CONNECTIVITY SECTION ID
101  LOGICAL :: READ_SEC1 ! FLAG FOR READING SECTION 1
102  LOGICAL :: READ_SEC2 ! FLAG FOR READING SECTION 2
103  LOGICAL :: READ_SEC3 ! FLAG FOR READING SECTION 3
104  INTEGER :: NELEMTOTAL ! TOTAL NUMBER OF UNV ELEMENTS
105  INTEGER :: NPOINT ! TOTAL NUMBER OF NODES
106  INTEGER :: NBFAMILY ! TOTAL NUMBER OF FAMILY
107  INTEGER :: NELIN ! TOTAL NUMBER OF INNER TRIANGLES
108  INTEGER :: SIZE_FLUX ! TOTAL NUMBER OF INNER SURFACES
109  INTEGER, DIMENSION(:), ALLOCATABLE :: VECTNB ! VECTEUR AUX POUR NACHB
110 
111  ! BEGIN MODIF V STOBIAC
112  ! FORMAT MED
113  INTEGER, PARAMETER :: MDIM=3 ! MESH DIMENSION (EXPECTED WITH ESTEL)
114 
115  CHARACTER(LEN=16), DIMENSION(MDIM) :: NAMECOOR, UNITCOOR ! NOM ET UNITE DES COORDONNES
116  CHARACTER(LEN=16) :: DATA_TEMP
117  CHARACTER(LEN=43) :: TXT, TXT_OLD
118 #if defined HAVE_MED
119  CHARACTER(LEN=MED_NAME_SIZE),DIMENSION(:),ALLOCATABLE :: FAM
120  CHARACTER(LEN=MED_NAME_SIZE) :: MED_VERSION ! VERSION OF THE MED FILE
121  CHARACTER(LEN=MED_NAME_SIZE) :: MESH_NAME ! CURRENT MESH NAME
122  CHARACTER(LEN=MED_NAME_SIZE) :: DTUNIT ! TIME STEP UNIT
123  CHARACTER(LEN=MED_COMMENT_SIZE) :: DESC ! MESH DESCRIPTION
124  CHARACTER(LEN=MED_LNAME_SIZE),DIMENSION(:),ALLOCATABLE ::
125  & gr_family
126 #endif
127  CHARACTER, DIMENSION(:), ALLOCATABLE :: DUMNAME
128 
129  LOGICAL :: HDFOK, MEDOK ! CHECK COMPATIBILITY
130 
131  INTEGER :: CRET ! ERROR CODE
132  INTEGER :: NBTRI2 ! TOTAL NB OF TRIANGLE(INNER AND BOUNDARY)
133  INTEGER :: NDIM ! DIMENSION OF THE PB (3)
134  INTEGER :: MAJOR, MINOR, REL ! MED VERSION OF THE FILE
135  INTEGER :: NBMESH ! TOTAL NUMBER OF MESH
136  INTEGER :: IDUM, NUM ! DUMMY
137  INTEGER :: STYPE ! MESH TYPE (CARTESIAN, UNSTRUCTURED)
138  INTEGER :: NBGRF
139  INTEGER :: FID ! FILE DESCRIPTOR FOR THE MED FILE
140  INTEGER :: NBFAMILY2 ! TEMPORARY NUMBER OF FAMILY
141  INTEGER, DIMENSION(:), ALLOCATABLE :: DUMNUM
142  INTEGER, DIMENSION(:), ALLOCATABLE :: IKLESTRI2
143  INTEGER, DIMENSION(:), ALLOCATABLE :: NUFATRIA, NUFATRIA2 ! TRIANGLE COLOR
144  INTEGER, DIMENSION(:), ALLOCATABLE :: NUFATETRA ! TETRAHEDRA COLOR
145  INTEGER, DIMENSION(:), ALLOCATABLE :: NUFANO ! NODES COLOR
146  INTEGER, DIMENSION(:,:), ALLOCATABLE :: ID_CHANGE_LOG
147  INTEGER, DIMENSION(:,:), ALLOCATABLE :: ID_CHANGE_LOG2
148 
149  DOUBLE PRECISION, ALLOCATABLE :: COOR(:) ! COORD NODES
150 
151  INTEGER :: PTET1,PTET2,PTET3,DEB1,FIN1,DEB2,FIN2,DEB3,FIN3
152  INTEGER :: PTRI1,PTRI2,PTRI3
153  LOGICAL :: FOUND_TET
154 
155  ! END MODIF V STOBIAC
156 
157  DOUBLE PRECISION, ALLOCATABLE :: X1(:),Y1(:),Z1(:) ! COORD NODES
158  INTEGER, ALLOCATABLE :: ECOLOR(:) ! ELEMENTS' COLOUR
159  INTEGER :: ELEM ! TYPE OF THE ELEMENT
160  INTEGER :: IKLE1,IKLE2,IKLE3,IKLE4,IKLEB ! NODES
161  INTEGER, DIMENSION(:), ALLOCATABLE :: IKLESTET ! CONNECTIVITE EN
162  ! RENUMEROTATION GLOBAL DE LA BIEF POUR LES TETRAEDRES
163  INTEGER, DIMENSION(:), ALLOCATABLE :: IKLESTRI ! CONNECTIVITE EN
164  ! RENUMEROTATION GLOBAL DE LA BIEF POUR LES TRIANGLES
165  INTEGER, DIMENSION(:,:), ALLOCATABLE :: IKLESTRIN ! CONNECTIVITE EN
166  ! RENUMEROTATION GLOBAL DE LA BIEF POUR LES TRIANGLES
167  INTEGER, DIMENSION(:,:), ALLOCATABLE :: IKLEIN ! COPIE AJUSTEE DE IKLESTRIN
168  INTEGER, DIMENSION(:,:), ALLOCATABLE :: TYPELEM ! TYPE D'ELT
169  INTEGER :: NBTET,NBTRI ! NBRE DE TETRA, TRIANGLE BORD
170  INTEGER, DIMENSION(:), ALLOCATABLE :: TETTRI, TETTRI2 ! JOINTURE
171  ! TETRA/TRIANGLE DE BORD
172  INTEGER, DIMENSION(:), ALLOCATABLE :: EPART ! NUMERO DE PARTITION
173  ! PAR ELEMENT
174  INTEGER, DIMENSION(:), ALLOCATABLE :: NPART ! NUMERO DE PARTITION
175  ! PAR NOEUDS
176  INTEGER, DIMENSION(:), ALLOCATABLE :: CONVTRI ! CONVERTISSEUR
177  ! NUMERO LOCAL TRIA/TETRA NUMERO GLOBAL; INVERSE DE TYPELEM(:,2)
178  INTEGER :: PARSEC ! RUNTIME
179  INTEGER, DIMENSION(:), ALLOCATABLE :: NPOINTSD, NELEMSD ! NBRE
180  ! DE POINTS ET D'ELEMENTS PAR SOUS-DOMAINE
181  INTEGER, DIMENSION(:), ALLOCATABLE :: NPOINTISD ! NBRE
182  ! DE POINTS D'INTERFACE PAR SOUS-DOMAINE
183  ! VECTEURS LIES AUX CONNECTIVITEES NODALES INVERSES
184  INTEGER, DIMENSION(:), ALLOCATABLE :: NODES1,NODES2,NODES3,NODES4
185  INTEGER, DIMENSION(:), ALLOCATABLE :: NODES1T,NODES2T,NODES3T
186  INTEGER, DIMENSION(:), ALLOCATABLE :: TRIUNV ! BUFFER POUR ECRIRE
187  ! DANS LES .UNV, D'ABORD LES TETRAS PUIS LES TRIA
188 ! POUR TRAITEMENT DES DIRICHLETS CONFONDUS AVEC L'INTERFACE
189  INTEGER :: NBCOLOR ! NBRE DE COULEUR DE MAILLES EXTERNES
190  INTEGER, DIMENSION(:), ALLOCATABLE :: PRIORITY
191  INTEGER, DIMENSION(:), ALLOCATABLE :: NCOLOR,NCOLOR2
192 ! POUR TRAITEMENT DES DIRICHLETS SUR LES NOEUDS DE TETRA
193  LOGICAL, DIMENSION(:,:), ALLOCATABLE :: TETCOLOR
194  LOGICAL, DIMENSION(:), ALLOCATABLE :: DEJA_TROUVE
195 ! INDISPENSABLE POUR PARALLELISME TELEMAC
196  INTEGER, DIMENSION(:), ALLOCATABLE :: KNOLG
197  INTEGER, DIMENSION(:,:), ALLOCATABLE :: NACHB
198  LOGICAL :: NACHBLOG
199 ! MAXIMUM GEOMETRICAL MULTIPLICITY OF A NODE (VARIABLE AUSSI
200 ! PRESENTE DANS LA BIEF, NE PAS CHANGER L'UNE SANS L'AUTRE)
201  INTEGER, PARAMETER :: NBMAXNSHARE = 10
202 ! CETTE VARIABLE EST LIEE A LA PRECEDENTE ET DIMENSIONNE DIFFERENTS
203 ! VECTEURS
204 ! NOTE SIZE OF NACHB WILL BE HERE 2 MORE THAN IN BIEF, BUT THE EXTRA 2 ARE
205 ! LOCAL WORK ARRAYS
206  INTEGER :: NBSDOMVOIS = nbmaxnshare + 2
207 !
208  INTEGER, PARAMETER :: MAX_SIZE_FLUX = 99
209 ! NUMBER OF INNER SURFACE (SAME AS SIZE_FLUX AT THE END)
210  INTEGER, DIMENSION(MAX_SIZE_FLUX) :: SIZE_FLUXIN
211 ! VECTEUR POUR PROFILING
212  INTEGER TEMPS_SC(20)
213 !
214 !F.D
215  INTEGER, DIMENSION(: ), ALLOCATABLE :: TEMPO,GLOB_2_LOC
216  INTEGER, DIMENSION(:,:), ALLOCATABLE :: IKLES,IKLE,IFABOR
217  INTEGER, DIMENSION(:,:), ALLOCATABLE :: NULONE,IKLBOR
218  INTEGER :: IKL,NSOLS,NSOLS_OLD,N1,N2
219  INTEGER :: IELEM,IPTFR,IELEB,N3
220  LOGICAL, DIMENSION(:), ALLOCATABLE :: FACE_CHECK
221  INTEGER, PARAMETER :: NCOL = 256
222  INTEGER, DIMENSION(NCOL ) :: COLOR_PRIO
223  INTEGER :: PRIO_NEW,NPTFR
224  INTEGER, DIMENSION(:), ALLOCATABLE :: NBOR2,NBOR
225  INTEGER, DIMENSION(:), ALLOCATABLE :: NELBOR,LIHBOR
226 
227 !D****************************************************** ADDED BY CHRISTOPHE DENIS
228  INTEGER, DIMENSION(:), ALLOCATABLE :: NELEM_P
229 ! SIZE NPARTS, NELEM_P(I) IS THE NUMBER OF FINITE ELEMENTS ASSIGNED TO SUBDOMAIN I
230  INTEGER, DIMENSION(:), ALLOCATABLE :: NPOIN_P
231 ! SIZE NPARTS, NPOIN_P(I) IS THE NUMBER OF NODES ASSIGNED TO SUBDOMAIN I
232  INTEGER :: NODE
233 ! ONE NODE ...
234  INTEGER :: POS_NODE
235 ! POSITION OF ONE ONE NODE
236  INTEGER :: MAX_NELEM_P
237 ! MAXIMUM NUMBER OF FINITE ELEMENTS ASSIGNED AMONG SUBDOMAINS
238  INTEGER :: MAX_NPOIN_P
239 ! MAXIMUM NUMBER OF NODES ASSIGNED AMONG SUBDOMAINS
240  INTEGER :: MAX_TRIA
241 ! MAXIMUM NUMBER OF TRIANGLE SHARING A NODE
242  INTEGER :: THE_TRI
243 ! ONE TRIANGLE
244  INTEGER :: JJ
245 ! INDEX COUNTER
246  INTEGER, DIMENSION(:), ALLOCATABLE :: NUMBER_TRIA
247 ! MAXIMUM NUMBER OF TRIANGLE SHARING A SAME NODE
248  INTEGER, DIMENSION(:,:), ALLOCATABLE :: ELEGL
249 ! SIZE MAX_NELEM_P,NPARTS, ELEGL(J,I) IS THE GLOBAL NUMBER OF LOCAL FINITE ELEMENT J IN SUBDOMAIN I
250  INTEGER, DIMENSION(:,:), ALLOCATABLE :: NODEGL
251 ! SIZE MAX_NPOIN_P,NPARTS, NODEGL(J,I) IS THE GLOBAL NUMBER OF LOCAL NODE J IN SUBDOMAIN I
252  INTEGER, DIMENSION(:,:), ALLOCATABLE :: NODELG
253 ! SIZE NPOINT, NODELG(J,I)=J, IS THE LOCAL NUMBER OF GLOBAL NODE J ON SUBDOMAIN I
254  INTEGER, DIMENSION(:,:), ALLOCATABLE :: TRI_REF
255 ! SIZE NPOINT*MAX_TRIA
256 ! EXTENS
257  CHARACTER(LEN=11) :: EXTENS
258  EXTERNAL extens
259 !D********************************************************
260  INTEGER SOMFAC(3,4)
261  parameter( somfac = reshape((/ 1,2,3 , 4,1,2 , 2,3,4 , 3,4,1 /),
262  & (/3,4/) ))
263 !
264 !-------------
265 ! 1. PREAMBULE
266 !-------------
267 !
268  CALL system_clock(count=temps_sc(1),count_rate=parsec)
269  ALLOCATE (vectnb(nbsdomvois-3))
270  WRITE(lu,*)' '
271  WRITE(lu,*)'+-------------------------------------------------+'
272  WRITE(lu,*)' PARTEL: TELEMAC ESTEL3D PARTITIONER'
273  WRITE(lu,*)'+-------------------------------------------------+'
274  WRITE(lu,*)' READING UNV AND LOG FILES'
275 !
276 ! WRITE FILE FORMAT AND TEST MED LIBRARY
277  IF (format_med) THEN
278 #if defined HAVE_MED
279  WRITE(lu,*) 'MED FORMAT IS USED FOR THE MESH'
280 #else
281  WRITE(lu,*) 'ERROR: MED FORMAT IS USED FOR THE MESH BUT MED
282  & LIBRARY IS NOT INSTALLED'
283  CALL plante(1)
284  stop
285 #endif
286  ELSE
287  WRITE(lu,*) 'UNV FORMAT IS USED FOR THE MESH'
288  ENDIF
289 
290 ! OPEN THE MESH FILE (MED)
291 #if defined HAVE_MED
292  IF (format_med) THEN
293 
294  ! CHECK IF THE FILE IS BOTH A MED & HDF5 FILE
295  CALL mficom (nameinp, hdfok, medok, cret)
296  CALL check_call(cret,'MFICOM')
297 
298  IF (.NOT.hdfok) WRITE(lu,*)'MESH FILE NOT COMPATIBLE WITH HDF5'
299  IF (.NOT.medok) WRITE(lu,*)'MESH FILE NOT COMPATIBLE WITH MED'
300 
301  ! OPEN
302  CALL mfiope(fid, nameinp, med_acc_rdonly, cret)
303  CALL check_call(cret,'MFIOPE')
304 
305  ! PRINT THE MED VERSION OF THE FILE
306  CALL mfisvr (fid, med_version, cret)
307  CALL check_call(cret,'MFISVR')
308  WRITE(lu,*) 'MED VERSION OF THE MESH FILE : '//trim(med_version)
309 
310  ! CHECK COMPATIBILITY (EXPECTED)
311  CALL mfinvr (fid, major, minor, rel, cret)
312  CALL check_call(cret,'MFINVR')
313  IF (major.LT.3) WRITE(lu,*) 'MED FILE IS TOO ''OLD'' (' //
314  & trim(med_version) // ') => PLEASE CONVERT IT WITH MEDIMPORT '
315  END IF
316 #endif
317 
318 ! OPEN THE MESH FILE (UNV)
319  IF (.NOT.format_med) THEN
320  OPEN(ninp,file=nameinp,status='OLD',form='FORMATTED',err=131)
321  rewind(ninp)
322  ENDIF
323 
324 ! OPEN THE COMPLEMENTARY FILE
325  OPEN(nlog,file=namelog,status='OLD',form='FORMATTED',err=130)
326  rewind(nlog)
327 ! END MODIF V STOBIAC
328 
329 !----------------------------------------------------------------------
330 ! 2A. LECTURE DU FICHIER .LOG
331 !---------------
332 ! The first line contains the number of nodes after a text descriptor.
333 ! We read a line, locate the colon ':' to then read the number.
334  READ(unit=nlog, fmt='(A200)', iostat=ios) line
335  IF (ios .NE. 0) THEN
336  WRITE(lu, *) 'ERROR READING THE MESH COMPLEMENTARY FILE.'
337  CALL plante(1)
338  ENDIF
339  pos =index(line,':') + 1
340  READ(unit=line(pos:),fmt=*, iostat=ios) npoint
341  IF (ios .NE. 0) THEN
342  WRITE(lu,*) 'FORMAT ERROR READING THE MESH COMPLEMENTARY FILE.'
343  CALL plante(1)
344  ENDIF
345 
346 ! The second line contains the number of elements after a text descriptor.
347 ! We read a line, locate the colon ':' and then read the number.
348  READ(unit=nlog, fmt='(A200)', iostat=ios) line
349  IF (ios .NE. 0) THEN
350  WRITE(lu,*)'ERROR READING THE MESH COMPLEMENTARY FILE.'
351  CALL plante(1)
352  ENDIF
353  pos =index(line,':') + 1
354  READ(unit=line(pos:),fmt=*, iostat=ios) nelemtotal
355  IF (ios .NE. 0) THEN
356  WRITE(lu,*)'FORMAT ERROR READING THE MESH COMPLEMENTARY FILE.'
357  CALL plante(1)
358  ENDIF
359 !
360 ! The third line contains the number of families after a text descripto.
361 ! We read a line, locate the colon ':' and then read the number.
362  READ(unit=nlog, fmt='(A200)', iostat=ios) line
363  IF (ios .NE. 0) THEN
364  WRITE(lu,*)'ERROR READING THE MESH COMPLEMENTARY FILE.'
365  CALL plante(1)
366  ENDIF
367  pos =index(line,':') + 1
368  READ(unit=line(pos:),fmt=*, iostat=ios) nbfamily
369  IF (ios .NE. 0) THEN
370  WRITE(lu,*)'FORMAT ERROR READING THE MESH COMPLEMENTARY FILE.'
371  CALL plante(1)
372  ENDIF
373 
374 ! BEGIN MODIF V STOBIAC
375 ! READ FAMILIES
376 #if defined HAVE_MED
377  IF (format_med) THEN
378  ! IGNORE NBFAMILY+1 LINES
379  DO i=1,nbfamily+1
380  READ(unit=nlog, fmt='(A200)', iostat=ios) line
381  IF (ios .NE. 0) THEN
382  WRITE(lu,*)'! PROBLEM WITH THE NUMBER OF FAMILY!'
383  ENDIF
384  ENDDO
385  END IF
386 #endif
387 
388  IF (.NOT.format_med) THEN
389  nbfamily=nbfamily+1 ! POUR TITRE DU BLOC
390  ALLOCATE(logfamily(nbfamily),stat=ierr)
391  CALL check_allocate(ierr,' LOGFAMILY')
392  DO i=1,nbfamily
393  READ(nlog,50,err=111,end=120)logfamily(i)
394  ENDDO
395  ENDIF
396 ! END MODIF V STOBIAC
397 
398  READ(unit=nlog, fmt='(A200)', iostat=ios) line
399  IF (ios .NE. 0) THEN
400  WRITE(lu,*)'! PROBLEM WITH THE NUMBER OF COLOR !'
401  CALL plante(1)
402  ENDIF
403  pos = index(line,':') + 1
404  READ(unit=line(pos:), fmt=*, iostat=ios) nbcolor
405  IF (ios .NE. 0) THEN
406  WRITE(lu,*)'! PROBLEM WITH THE NUMBER COLOR !'
407  ENDIF
408 
409  ALLOCATE(priority(nbcolor),stat=ierr)
410  CALL check_allocate(ierr,' PRIORITY')
411  WRITE(lu,92) npoint
412  WRITE(lu,93) nelemtotal
413  WRITE(lu,94) nbcolor
414  IF (nbcolor.EQ.0) THEN
415  WRITE(lu,*) 'VOUS AVEZ OUBLIE DE REMPLIR LE FICHIER LOG...'
416  CALL plante(1)
417  stop
418  ENDIF
419 
420  ! MODIFICATION JP RENAUD 15/02/2007
421  ! SOME TEXT HAS BEEN ADDED BEFORE THE LIOST OF PRIORITIES.
422  ! READ A 200 CHARACTER LINE, FIND THE ':' AND THEN
423  ! READ THE VALUES AFTER THE ':'
424  READ(unit=nlog, fmt='(A200)', iostat=ios) line
425  IF (ios .NE. 0) THEN
426  ! '!------------------------------------------!'
427  WRITE(lu,*)'! PROBLEM WITH THE PRIORITY OF COLOR NODES !'
428  CALL plante(1)
429  ENDIF
430 
431  pos = index(line,':') + 1
432  READ(unit=line(pos:), fmt=*, iostat=ios) (priority(j),j=1,nbcolor)
433  IF (ios .NE. 0) THEN
434  ! '!------------------------------------------!'
435  WRITE(lu,*)'! PROBLEM WITH THE PRIORITY OF COLOR NODES !'
436  CALL plante(1)
437  ENDIF
438  ! END MODIFICATION JP RENAUD
439  WRITE(lu,*) (priority(j),j=1,nbcolor)
440  CLOSE(nlog)
441 !
442 ! 2B. ALLOCATIONS MEMOIRES ASSOCIEES
443 !---------------
444 
445 !D ****************************** ALLOCATION MEMORY ADDED BY CD
446  ALLOCATE(nelem_p(nparts),stat=ierr)
447  CALL check_allocate(ierr,'NELEM_P')
448  ALLOCATE(npoin_p(nparts),stat=ierr)
449  CALL check_allocate(ierr,'NPOIN_P')
450 !D *******************************
451 
452 ! BEGIN MODIF VSTOBIAC
453 ! READ THE MESH FILE (MED)
454 #if defined HAVE_MED
455  IF (format_med) THEN
456 
457  ! DEBUG MODE : FOR DEBUGGING CRASH IN MED ROUTINES
458  ! READING THE NUMBER OF MESH IN THE FILE
459  CALL mmhnmh(fid, nbmesh, cret)
460  CALL check_call(cret,'MMHNMH')
461  IF (nbmesh .NE. 1) WRITE(lu,*)'! ONLY ONE MESH EXPECTED !'
462 
463  ! LECTURE DES INFOS SUR LE MAILLAGE
464  CALL mmhmii(fid,1,mesh_name,ndim,idum,stype,desc,dtunit,idum,
465  & idum,idum,namecoor,unitcoor,cret)
466  CALL check_call(cret,'MMHMII')
467 
468  ! GET NUMBER OF FAMILIES IN THE MESH
469  CALL mfanfa(fid,mesh_name,nbfamily2,cret)
470  CALL check_call(cret,'MFANFA')
471 
472  ! AND THEN ALLOCATE...
473  ALLOCATE (fam(nbfamily2))
474 
475  ! READING FAMILIES, GROUPS, ATTRIBUTES IN THE MESH
476  ALLOCATE(id_change_log2(nbfamily2,2))
477  id_change_log2 = 0
478  nbfamily = 0
479 
480  DO i=1,nbfamily2
481  ! THERE CAN BE MORE THAN ONE GROUP
482  ! WE CONSIDER THERE ARE AT MAX 10 GROUPS
483  !
484  ! GET THE NUMBER OF GROUPS IN THE FAMILIY
485  CALL mfanfg (fid,mesh_name,i,nbgrf,cret)
486  CALL check_call(cret,'MFANFG')
487 
488  IF (nbgrf.NE.0) THEN
489  !
490  ! AND THEN ALLOCATE...
491  ALLOCATE(gr_family(nbgrf))
492  !
493  ! /! MED 3 ONLY !\ : GETTING INFORMATION FOR THE FAMILY
494  CALL mfafai(fid,mesh_name,i,fam(i),num,gr_family,
495  & cret)
496  CALL check_call(cret,'MFAFAI')
497  !
498  ! WARNING (ASSUMPTION) : WE EXPECT ONLY 1 GROUP PER FAMILY...
499  ! WARNING (ASSUMPTION) : SYNTAX EXPECTED : <COLOR_ID>:<NAME_GROUP>
500  pos = index(gr_family(1),':')-1
501  !
502  ! BEGIN MODIF V STOBIAC
503  IF ((pos.LT.1).OR.(pos.GT.3)) THEN
504  ! COULD BE AN ERROR OR A GROUP OF NODE
505  ! (MAY HAPPEN WHEN THE MED FILE IS EXPORTED FROM AN UNV FILE)
506  pos = 0
507  nsols = -99
508 
509 ! WRITE(LU,*)'WARNING : ERROR WHEN READING GROUP NUMBER'
510 ! & // GR_FAMILY(1)//' SYNTAX : <COLOR_ID>:'
511 ! & // '<NAME_GROUP>'
512 
513  ELSE
514  !
515  ! THE FORMAT IS ADAPTED TO POS IN ORDER TO
516  ! PREVENT PROBLEMS ON THE BGQ CLUSTER
517  READ(gr_family(1)(1:pos),fmt='(i'//achar(48+pos)//')',
518  & iostat=ios) nsols
519  ! END MODIF V STOBIAC
520  IF (ios.NE.0) THEN
521  nsols = -99
522 ! WRITE(LU,*)'CONVERTING ERROR, FROM GROUP:'//
523 ! & TRIM(GR_FAMILY(1))//'TO COLOR:',NSOLS
524  ENDIF
525  ENDIF
526  !
527  ! TABLE FOR INDEXING MED_ID (INTERN TO MED) TO USER'S COLOR
528  id_change_log2(i,1) = nsols
529  id_change_log2(i,2) = num
530  !
531  ! COUNT THE NUMBER OF USEFULL FAMILY
532  IF ((nsols.NE.-99).AND.(nsols.GT.0))
533  & nbfamily = nbfamily + 1
534  DEALLOCATE(gr_family)
535  ENDIF
536  ENDDO
537 
538  ! PRINT A SUMMARY
539  print*,'NBFAMILY2',nbfamily2
540  DO i=1,nbfamily2
541  WRITE(lu,*) 'MED_ID=',id_change_log2(i,2),
542  & ' <==> NCOLOR_ID=',id_change_log2(i,1)
543  ENDDO
544 
545  ! CORRECTION OF ID_CHANGE_LOG TO ERASE USELESS FAMILY
546  ALLOCATE(id_change_log(nbfamily,2))
547  j = 0
548  DO i = 1, nbfamily2
549  nsols = id_change_log2(i,1)
550  IF ((nsols.NE.-99).AND.(nsols.GT.0)) THEN
551  j = j + 1
552  id_change_log(j,:) = id_change_log2(i,:)
553  ENDIF
554  ENDDO
555  DEALLOCATE(id_change_log2)
556  DEALLOCATE(fam)
557 
558  ! PRINT A SUMMARY
559  print*,'NBFAMILY',nbfamily
560  DO i=1,nbfamily
561  WRITE(lu,*) 'MED_ID=',id_change_log(i,2),' <==> NCOLOR_ID=',
562  & id_change_log(i,1)
563  ENDDO
564 
565  ! READ MESH ENTITIES (AXIS, NODES)
566  !-------------------------------------------------------
567 
568  ! READING THE NUMBER OF AXIS (HAVE BEEN ALREADY READEN IN MMHMII)
569  CALL mmhnax(fid,1,ndim,cret)
570  CALL check_call(cret,'MMHNAX')
571 
572  IF (ndim.NE.mdim) THEN ! E.G. = 3
573  WRITE(lu,*) '3 DIMENSIONS MESH EXPECTED - FOUND ', ndim
574  END IF
575 
576  ! READING THE NUMBER OF NODES
577  CALL mmhnme (fid, mesh_name, med_no_dt, med_no_it, med_node,
578  & med_none, med_coordinate, med_nodal, idum, idum,
579  & npoint, cret)
580  CALL check_call(cret,'MMHNME')
581 
582  ! ALLOCATING BEFORE READING THE COORDINATES
583  ALLOCATE(coor(npoint*mdim), stat=ierr)
584  CALL check_allocate(ierr,' COOR')
585 
586  CALL mmhcor(fid,trim(mesh_name),med_no_dt,med_no_it,
587  & med_no_interlace, coor, cret)
588  CALL check_call(cret,'MMHCOR')
589 
590  ALLOCATE(x1(npoint), stat=ierr)
591  CALL check_allocate(ierr,' X1')
592  ALLOCATE(y1(npoint), stat=ierr)
593  CALL check_allocate(ierr,' Y1')
594  ALLOCATE(z1(npoint), stat=ierr)
595  CALL check_allocate(ierr,' Z1')
596 
597  DO i = 1, npoint
598  x1(i) = coor(i)
599  y1(i) = coor(i+ npoint)
600  z1(i) = coor(i+2*npoint)
601  END DO
602  DEALLOCATE(coor)
603 
604  ! READ MESH ENTITIES (2D & 3D ELEMENTS)
605  !---------------------------------------
606 
607  ! ATTENTION : DIFFERENT FROM UNV FORMAT
608  ! GET TOTAL NUMBER OF TRIANGLES (TOTAL = INNER + BORDER MESH)
609  CALL mmhnme(fid,mesh_name,med_no_dt,med_no_it,med_cell,
610  & med_tria3, med_connectivity,med_nodal,idum,idum,
611  & nbtri2,cret)
612  CALL check_call(cret,'MMHNME')
613 
614  ! GET NUMBER OF TETRA (3D MESH)
615  CALL mmhnme(fid,mesh_name,med_no_dt,med_no_it,med_cell,
616  & med_tetra4, med_connectivity,med_nodal,idum,idum,
617  & nbtet,cret)
618  CALL check_call(cret,'MMHNME')
619 
620  ! ATTENTION : NBTRI = INNER + BORDER TRIANGLE
621  ! PRINT A SUMMARY
622  WRITE(lu,*)'NUMBER OF TETRAHEDRONS IN THE MESH',nbtet
623  WRITE(lu,*)'NUMBER OF TRIANGLES IN THE MESH', nbtri2
624  WRITE(lu,*)'NUMBER OF POINTS IN THE MESH', npoint
625  WRITE(lu,*)'NUMBER OF EXTERNAL FACES', nbcolor
626 
627  txt = ' '
628  DO i=1, nbcolor
629  data_temp = i2char2(priority(i))
630  txt_old = txt
631  txt = trim(txt_old) // ' ' // trim(data_temp)
632  ENDDO
633 
634  WRITE(lu,*)'PRIORITE DES FACES EXTERNES'//trim(txt)
635 
636  ALLOCATE(iklestet(nbtet*4), stat=ierr)
637  CALL check_allocate(ierr,' IKLESTET')
638  ALLOCATE(nufatetra(nbtet), stat=ierr)
639  CALL check_allocate(ierr,' NUFATETRA')
640  ALLOCATE(dumname(nbtet),stat=ierr)
641  CALL check_allocate(ierr,' DUMNAME')
642  ALLOCATE(dumnum(nbtet),stat=ierr)
643  CALL check_allocate(ierr,' DUMNUM')
644 
645  ! READ TETRA CONNECTIVITY AND TETRA FAMILIY
646  CALL mmhelr(fid,mesh_name,med_no_dt,med_no_it,med_cell,
647  & med_tetra4,med_nodal,med_no_interlace,iklestet,idum,
648  & dumname,idum,dumnum,idum,nufatetra,cret)
649  CALL check_call(cret,'MMHELR')
650  DEALLOCATE(dumname,dumnum)
651 
652  ! YA_PARA_MODIF
653  DO i=1,nbtet
654  DO j=1,nbfamily
655  IF (nufatetra(i).EQ.id_change_log(j,2)) THEN
656  nufatetra(i) = id_change_log(j,1)
657  EXIT
658  ENDIF
659  ENDDO
660  ENDDO
661  ! YA_PARA_MODIF FIN
662 
663  ! READ TRIANGLES CONNECTIVITY AND TRIANGLES FAMILIY
664  ! S. FALAPPI : NBTRI2 MAY BE NULL IF WE TRY TO READ A MESH FILE
665  ! WRITTEN BY ESTEL-3D WHERE NO BORDER CONNECTIVITY IS WRITTEN (DEBUG PURPOSES)
666  IF (nbtri2.GT.0) THEN
667 
668  ALLOCATE(iklestri2(nbtri2*3), stat=ierr)
669  CALL check_allocate(ierr,' IKLESTRI2')
670  ALLOCATE(nufatria2(nbtri2), stat=ierr)
671  CALL check_allocate(ierr,' NUFATRIA2')
672  ALLOCATE(dumname(nbtri2),stat=ierr)
673  CALL check_allocate(ierr,' DUMNAME')
674  ALLOCATE(dumnum(nbtri2),stat=ierr)
675  CALL check_allocate(ierr,' DUMNUM')
676 
677  CALL mmhelr(fid,mesh_name,med_no_dt,med_no_it,med_cell,
678  & med_tria3,med_nodal,med_no_interlace,iklestri2,idum,
679  & dumname,idum,dumnum,idum,nufatria2,cret)
680  CALL check_call(cret,'MMHELR')
681  DEALLOCATE(dumname,dumnum)
682 
683  ALLOCATE(iklestrin(nbtri2,3), stat=ierr)
684  CALL check_allocate(ierr,' IKLESTRIN')
685 
686  ! CORRECTION OF ORDER TO MATCH UNV FORMAT
687  DO i = 1, nbtri2
688  iklestrin(i,1) = iklestri2(i)
689  iklestrin(i,2) = iklestri2(i+nbtri2)
690  iklestrin(i,3) = iklestri2(i+2*nbtri2)
691  ENDDO
692  DO i = 1, nbtri2
693  iklestri2(3*(i-1)+1:3*i) = iklestrin(i,1:3)
694  ENDDO
695  DEALLOCATE(iklestrin)
696 
697  DO i=1,nbtri2
698  DO j=1,nbfamily
699  IF (nufatria2(i).EQ.id_change_log(j,2)) THEN
700  nufatria2(i) = id_change_log(j,1)
701  EXIT
702  ENDIF
703  ENDDO
704  ENDDO
705  ENDIF
706 
707  !-----------------------------------------------------------------------
708  ! BUILD NUFANO (BORDER MESH COLOR)
709  ! On va reconstruire nufano a partir de la connectivite et des priorites.
710  ! On part de la priorite la moins elevee (la derniere) jusqu'a la premiere
711  ! Petit test pour savoir si on relit un maillage med de sortie de code ou
712  ! un maillage med cree de toute piece.
713  ! Dans le premier cas, on a aucune information sur les faces
714  ! et par contre nufano est juste
715  ALLOCATE(nufano(npoint), stat=ierr)
716  CALL check_allocate(ierr,' NUFANO')
717 
718  IF (nbtri2.GT.0) THEN
719  nufano = 0
720  DO i=nbcolor,1,-1
721  DO j=1,nbtri2
722  IF (nufatria2(j)==priority(i)) THEN
723  nufano(iklestri2(3*(j-1)+1)) = priority(i)
724  nufano(iklestri2(3*(j-1)+2)) = priority(i)
725  nufano(iklestri2(3*j)) = priority(i)
726  ENDIF
727  ENDDO
728  ENDDO
729 
730  size_fluxin(:) = 0
731  size_flux = 0
732  nelin = 0
733  nbtri = 0
734 
735  DO j=1,nbtri2
736 
737  IF (nufatria2(j).GT.0) THEN
738 
739  ! NSOLS_OLD IS USED FOR SAVING USE OF A NEW VARIABLE
740  nsols_old = nufatria2(j)
741 
742  ! INNER SURFACE ARE SUPPOSED TO HAVE A NSOLS VALUE > 100
743  IF (nsols_old.GT.maxval(priority)) THEN
744  nsols_old = nsols_old - 100 + maxval(priority)
745  ENDIF
746 
747  prio_new = size_fluxin(nsols_old)
748  IF (prio_new.EQ.0) THEN
749  size_flux = size_flux + 1
750  size_fluxin(nsols_old) = 1
751  ENDIF
752 
753  IF (nufatria2(j).LT.100) nbtri = nbtri + 1
754 
755  IF (nufatria2(j).GE.100) nelin = nelin + 1
756  ENDIF
757  ENDDO
758 
759  ALLOCATE(iklestrin(nelin,4))
760  CALL check_allocate(ierr,' IKLESTRIN')
761  iklestrin(:,:)=-999
762 
763  i = 0
764  DO j=1,nbtri2
765 
766  IF (nufatria2(j).GE.100) THEN
767  i = i + 1
768  ikle1 = iklestri2(3*(j-1)+1)
769  ikle2 = iklestri2(3*(j-1)+2)
770  ikle3 = iklestri2(3*j)
771 
772  iklestrin(i,1) = nufatria2(j)
773  iklestrin(i,2) = ikle1
774  iklestrin(i,3) = ikle2
775  iklestrin(i,4) = ikle3
776  ENDIF
777  ENDDO
778  ENDIF
779 
780  ! WRITE THE REAL NUMBER OF BORDER AND INNER TRIANGLES
781  print*,'CORRECTED NB OF BORDER TRIANGLES NBTRI',nbtri
782  print*,'NUMBER OF INNER TRIANGLES NELIN',nelin
783 
784  ! CORRECTION OF THE TOTAL NUMBER OF ELEMENT
785  nelemtotal = nbtet + nbtri
786  print*,'CORRECTED TOTAL NUMBER OF ELEMENT',nelemtotal
787 
788  ! ALLOCTAION OF IKLESTRI WHICH MATCH UNV FORMAT
789  ALLOCATE(iklestri(3*nbtri))
790  CALL check_allocate(ierr,' IKLESTRI')
791  ALLOCATE(nufatria(nbtri), stat=ierr)
792  CALL check_allocate(ierr,' NUFATRIA')
793 
794  ! CORRECTION OF IKLESTRI TO REMOVE INNER TRIANGLES
795  i = 0
796  DO j=1,nbtri2
797  ! ONLY BORDER TRIANGLES
798  IF ((nufatria2(j).GT.0) .AND. (nufatria2(j).LT.100)) THEN
799  i = i + 1
800  iklestri(3*(i-1)+1) = iklestri2(3*(j-1)+1)
801  iklestri(3*(i-1)+2) = iklestri2(3*(j-1)+2)
802  iklestri(3*i) = iklestri2(3*j)
803  nufatria(i) = nufatria2(j)
804  ENDIF
805  ENDDO
806  DEALLOCATE(iklestri2)
807  DEALLOCATE(nufatria2)
808 
809  ALLOCATE(ikle(nbtet,4), stat=ierr)
810  CALL check_allocate(ierr,' IKLE')
811  ikle = 0
812  DO i = 1, nbtet
813  ! CORRECTION OF ORDER TO MATCH UNV FORMAT
814  ikle(i,1) = iklestet(i+nbtet)
815  ikle(i,2) = iklestet(i)
816  ikle(i,3) = iklestet(i+2*nbtet)
817  ikle(i,4) = iklestet(i+3*nbtet)
818  ENDDO
819 
820  ALLOCATE(typelem(nelemtotal,2),stat=ierr)
821  CALL check_allocate(ierr,' TYPELEM')
822 
823  ! RESCAN TO CHANGE ORDER IN IKLESTET
824  ! AND FILE ELEMENT TABLE
825  DO i = 1, nbtet
826  iklestet(4*(i-1)+1:4*i) = ikle(i,1:4)
827  typelem(i,1) = 111
828  typelem(i,2) = i
829  ENDDO
830 
831  ALLOCATE(nbor2(npoint), stat=ierr)
832  CALL check_allocate(ierr,' NBOR2')
833  ALLOCATE(glob_2_loc(npoint))
834  CALL check_allocate(ierr,' GLOB_2_LOC')
835  glob_2_loc(:) = 0 ! GLOBAL TO LOCAL NUMBERING
836  iptfr = 0
837 
838  ! DEFINITION OF BOUNDARY NODE COLOR,
839  ! GLOBAL TO LOCAL CONVERTER AND
840  ! ELEMENT TABLE
841  DO i = 1, nbtri
842  DO j = 1, 3
843  ikl = iklestri(3*(i-1)+j)
844  k = glob_2_loc(ikl)
845  IF ((k.EQ.0).AND.(nufano(ikl)>0)) THEN
846  iptfr = iptfr + 1
847  nbor2(iptfr) = ikl
848  glob_2_loc(ikl) = iptfr
849  ENDIF
850  ENDDO
851  ENDDO
852  nptfr = iptfr
853  DEALLOCATE(glob_2_loc)
854 
855  ALLOCATE(convtri(nelemtotal), stat=ierr)
856  CALL check_allocate(ierr,' CONVTRI')
857  ! DEFINITION OF TYPE FOR TRIANGLE
858  DO i = 1,nbtri
859  idum = nbtet + i
860  convtri(i) = idum
861  typelem(idum,1) = 91
862  typelem(idum,2) = i
863  ENDDO
864 
865  ALLOCATE(ecolor(nelemtotal), stat=ierr)
866  CALL check_allocate(ierr,' ECOLOR')
867 
868  ! DEFINITION OF ELEMENT COLOR (TETRA AND TRIANGLE)
869  ecolor = 0
870  DO i = 1, nbtet
871  ecolor(i) = nufatetra(i)
872  END DO
873  DEALLOCATE(nufatetra)
874  DO i = 1, nbtri
875  ecolor(i+nbtet) = nufatria(i)
876  END DO
877  DEALLOCATE(nufatria)
878 
879  ALLOCATE(ncolor2(npoint), stat=ierr)
880  CALL check_allocate(ierr,' NCOLOR2')
881 
882  ncolor2(:) = 0
883  DO i = 1, nptfr
884  ncolor2(nbor2(i)) = nufano(nbor2(i))
885  END DO
886  DEALLOCATE(nufano)
887 
888  ! ALLOCATION REQUIRED FOR THE FUTUR
889  ALLOCATE(npointsd(nparts),stat=ierr)
890  CALL check_allocate(ierr,' NPOINTSD')
891  ALLOCATE(nelemsd(nparts),stat=ierr)
892  CALL check_allocate(ierr,' NELEMSD')
893  ALLOCATE(npointisd(nparts),stat=ierr)
894  CALL check_allocate(ierr,' NPOINTISD')
895  ENDIF
896 
897 #endif
898 
899 ! READ THE MESH FILE (UNV)
900  IF (.NOT.format_med) THEN
901 
902  ! ALLOCATIONS
903  ALLOCATE(x1(npoint),stat=ierr)
904  CALL check_allocate(ierr,' X1')
905  ALLOCATE(y1(npoint),stat=ierr)
906  CALL check_allocate(ierr,' Y1')
907  ALLOCATE(z1(npoint),stat=ierr)
908  CALL check_allocate(ierr,' Z1')
909  ALLOCATE(ncolor(npoint),stat=ierr)
910  CALL check_allocate(ierr,' NCOLOR')
911  ALLOCATE(ncolor2(npoint),stat=ierr)
912  CALL check_allocate(ierr,' NCOLOR2')
913  ALLOCATE(ecolor(nelemtotal),stat=ierr)
914  CALL check_allocate(ierr,' ECOLOR')
915  ALLOCATE(iklestet(4*nelemtotal),stat=ierr)
916  CALL check_allocate(ierr,' IKLESTET')
917  ALLOCATE(iklestri(3*nelemtotal),stat=ierr)
918  CALL check_allocate(ierr,' IKLESTRI')
919  ALLOCATE(iklestrin(nelemtotal,4),stat=ierr)
920  CALL check_allocate(ierr,' IKLESTRIN')
921  ALLOCATE(typelem(nelemtotal,2),stat=ierr)
922  CALL check_allocate(ierr,' TYPELEM')
923  ALLOCATE(convtri(nelemtotal),stat=ierr)
924  CALL check_allocate(ierr,' CONVTRI')
925  ALLOCATE(npointsd(nparts),stat=ierr)
926  CALL check_allocate(ierr,' NPOINTSD')
927  ALLOCATE(nelemsd(nparts),stat=ierr)
928  CALL check_allocate(ierr,' NELEMSD')
929  ALLOCATE(npointisd(nparts),stat=ierr)
930  CALL check_allocate(ierr,' NPOINTISD')
931 
932 !F.D
933  ALLOCATE(nbor2(npoint),stat=ierr)
934  CALL check_allocate(ierr,' NBOR2')
935  ALLOCATE(tempo(2*npoint),stat=ierr)
936  CALL check_allocate(ierr,' TEMPO')
937  ALLOCATE(face_check(nbfamily),stat=ierr)
938  CALL check_allocate(ierr,' FACE_CHECK')
939  ALLOCATE(glob_2_loc(npoint),stat=ierr)
940  CALL check_allocate(ierr,' GLOB_2_LOC')
941  ALLOCATE(ikles(nelemtotal,4),stat=ierr)
942 
943  read_sec1 = .true.
944  read_sec2 = .true.
945  read_sec3 = .true.
946 
947  DO WHILE ( read_sec1 .OR. read_sec2 .OR. read_sec3 )
948 
949  moins1 = ' '
950  blanc = '1111'
951  DO WHILE (moins1/='-1' .OR. blanc/=' ')
952  READ(ninp,2000, err=1100, end=1200) blanc, moins1
953  END DO
954 
955  2000 FORMAT(a4,a2)
956 
957  nsec = -1
958 
959  DO WHILE (nsec == -1)
960  READ(ninp,*, err=1100, end=1200) nsec
961  END DO
962 
963  SELECT CASE (nsec)
964 
965  CASE ( nsec1 )
966 
967  read_sec1 = .false.
968 
969  READ(ninp,25,err=1100, end=1200) titre
970 
971  25 FORMAT(a80)
972 
973  CASE ( nsec2 )
974 
975  read_sec2 = .false.
976  ncolor(:) = -1
977  tempo(:) = 0
978 
979  DO ielem = 1, npoint
980  READ(ninp,*,err=1100,end=1200) n,n1,n2,ncolor(ielem)
981  READ(ninp,*,err=1100,end=1200) x1(ielem), y1(ielem),
982  & z1(ielem)
983  tempo(n) = ielem
984  ENDDO
985 
986  CASE (nsec3 )
987 
988  read_sec3 = .false.
989 
990  nbtet = 0 ! NUMBER OF TETRA ELEMENTS TO 0
991  nbtri = 0 ! NUMBER OF BORDER ELEMENTS TO 0
992  nptfr = 0 ! NUMBER OF BORDER NODES TO 0.
993  nelin = 0 ! NUMBER OF INNER SURFACES TO 0.
994  size_flux = 0 ! NUMBER OF USER SURFACES TO 0.
995  nbor2(:) = 0 ! LOCAL TO GLOBAL NUMBERING
996  glob_2_loc(:) = 0 ! GLOBAL TO LOCAL NUMBERING
997 
998 !OB'S STFF
999  ecolor(:) = -1
1000  iklestet(:) = -1
1001  iklestri(:) = -1
1002  typelem(:,:) = -1
1003  convtri(:) = -1
1004 !
1005  iklestrin(:,:) = -1
1006 
1007  face_check(:) = .false.
1008  !
1009  color_prio(:) = 0
1010  size_fluxin(:) = 0
1011  !
1012  DO k = 1, nbcolor
1013  color_prio(priority(k)) = k
1014  END DO
1015 
1016  DO ielem = 1, nelemtotal
1017 
1018  READ(ninp,*,err=1100,end=1200) nsec, elem, n1, n2,
1019  & nsols,n3
1020 
1021  IF (nsec == -1) EXIT
1022 
1023  SELECT CASE ( elem )
1024 
1025  CASE ( 111 )
1026 
1027  nbtet = nbtet + 1
1028 
1029  ecolor(ielem) = nsols
1030 
1031  READ(ninp,*, err=1100, end=1200) ikle1, ikle2,
1032  & ikle3,ikle4
1033 
1034  ikles(ielem, 1) = tempo(ikle1)
1035  ikles(ielem, 2) = tempo(ikle2)
1036  ikles(ielem, 3) = tempo(ikle3)
1037  ikles(ielem, 4) = tempo(ikle4)
1038 
1039 !OB'S STFF
1040  n=4*(nbtet-1)+1
1041  iklestet(n)=ikle1 ! VECTEUR DE CONNECTIVITE
1042  iklestet(n+1)=ikle2
1043  iklestet(n+2)=ikle3
1044  iklestet(n+3)=ikle4
1045  typelem(ielem,1)=elem ! POUR TYPER LES ELTS
1046  typelem(ielem,2)=nbtet ! POUR CONVERSION NUM ELT> NUM TETRA
1047 
1048  CASE ( 91 )
1049 
1050  IF (nsols.GT.0.AND.nsols.LT.100) THEN
1051 
1052  IF ( nsols > ncol ) THEN
1053  WRITE(lu,*) 'COLOR ID POUR SURFACES EXTERNE'
1054  & // 'S TROP GRAND. LA LIMITE EST : ',ncol
1055  END IF
1056 
1057  prio_new = color_prio(nsols)
1058 
1059  IF ( prio_new .EQ. 0 ) THEN
1060  WRITE(lu,*) ' NUMERO DE FACE NON DECLARE',
1061  & 'DANS LE TABLEAU UTILISATEUR LOGFAMILY ',
1062  & 'VOIR LE FICHIER DES PARAMETRES '
1063  CALL plante(1)
1064  END IF
1065 
1066  face_check(prio_new) = .true.
1067 
1068  nbtri = nbtri + 1
1069 
1070  ecolor(ielem) = nsols
1071 
1072  READ(ninp,*,err=1100,end=1200)ikle1,ikle2,ikle3
1073  !
1074  prio_new = size_fluxin(nsols)
1075  !
1076  IF (prio_new.EQ.0) THEN
1077  size_flux = size_flux + 1
1078  size_fluxin(nsols) = 1
1079  ENDIF
1080 
1081  ikles(ielem, 1) = tempo(ikle1)
1082  ikles(ielem, 2) = tempo(ikle2)
1083  ikles(ielem, 3) = tempo(ikle3)
1084 
1085 !OB'S STFF
1086  n=3*(nbtri-1)+1
1087  iklestri(n)=ikle1
1088  iklestri(n+1)=ikle2
1089  iklestri(n+2)=ikle3
1090  typelem(ielem,1)=elem ! IDEM QUE POUR TETRA
1091  typelem(ielem,2)=nbtri
1092  convtri(nbtri)=ielem
1093 
1094  DO j=1,3
1095 
1096  ikl = ikles(ielem,j)
1097 
1098  iptfr = glob_2_loc(ikl)
1099 
1100  IF ( iptfr .EQ. 0 ) THEN
1101 
1102  nptfr = nptfr+1
1103  nbor2(nptfr) = ikl
1104  glob_2_loc(ikl) = nptfr
1105  iptfr = nptfr
1106 
1107  END IF
1108 
1109  ENDDO ! LOOP OVER THE NODES OF THE ELEMENT
1110 
1111  ELSE IF (nsols.GT.100) THEN
1112  !
1113  ! USER-DEFINED SURFACE FOR FLUXES COMPUTATION
1114  !
1115  ! NELIN IS THE COUNTER FOR THE INTERNAL ELEMENTS.
1116  ! ACTUALLY, WE ARE READING THE NEXT INTERNAL ELEMENT.
1117 
1118  ! NSOLS_OLD IS USED FOR SAVING USE OF A NEW VARIABLE
1119  nsols_old = nsols
1120  !
1121  ! PRIO_NEW IS USED FOR SAVING USE OF A NEW VARIABLE
1122  prio_new = size_fluxin(nsols_old)
1123  !
1124  IF (prio_new.EQ.0) THEN
1125  size_flux = size_flux + 1
1126  size_fluxin(nsols_old) = 1
1127  ENDIF
1128  !
1129  nelin = nelin + 1
1130  !
1131  READ(ninp,*,err=1100,end=1200)ikle1,ikle2,ikle3
1132  !
1133  iklestrin(nelin,1) = nsols
1134  iklestrin(nelin,2) = tempo(ikle1)
1135  iklestrin(nelin,3) = tempo(ikle2)
1136  iklestrin(nelin,4) = tempo(ikle3)
1137  !
1138  ELSE ! THIS IS AN INNER SURFACE, JUST READ THE LINE.
1139 
1140  READ(ninp,*,err=1100,end=1200)ikle1,ikle2,ikle3
1141 
1142  END IF
1143 
1144  CASE DEFAULT ! THIS IS AN UNKNOWN ELEMENT.
1145 
1146  WRITE(lu,*) 'ELEMENT INCONNU DANS LE MAILLAGE'
1147 
1148  END SELECT ! THE TYPE OF THE MESH ELEMENT
1149 
1150  END DO ! LOOP OVER ELEMENTS TO READ.
1151 
1152  DO k=1,nbcolor
1153  IF (.NOT. face_check(k)) THEN
1154  WRITE(lu,*) ' LA COULEUR DE FACE ',k,
1155  & ' N''APPARAIT PAS DANS LE MAILLAGE.'
1156  END IF
1157  END DO
1158 
1159 !-----------------------------------------------------------------------
1160 
1161  END SELECT ! TYPE OF THE SECTION
1162 
1163  END DO ! WHILE LOOP OVER SECTIONS TO READ
1164 
1165 !------------------------------------------------------- FIN VERSION F.D
1166 
1167  ENDIF
1168 
1169  ! CORRECTION DU NOMBRE D'ELEMENTS TOTAL CAR CELUI DANS LE .LOG EST
1170  ! COMPORTE DES ELEMENTS NON PRIS EN COMPTE DANS UNE ETUDE ESTEL
1171  nelemtotal=nbtet+nbtri
1172 
1173 ! FIN MODIF V STOBIAC
1174 
1175  CALL system_clock(count=temps_sc(2),count_rate=parsec)
1176  WRITE(lu,*)' TEMPS DE LECTURE FICHIERS LOG & UNV',
1177  & (1.0*(temps_sc(2)-temps_sc(1)))/(1.0*parsec),' SECONDS'
1178 
1179 !----------------------------------------------------------------------
1180 ! 3A. CONSTRUCTION DE TETTRI/TETTRI2: CORRESPONDANCE TETRA > TRIA
1181 !---------------
1182 
1183  ALLOCATE(nelbor(nbtri),stat=ierr)
1184  CALL check_allocate(ierr,' NELBOR')
1185  ALLOCATE(iklbor(nbtri,3),stat=ierr)
1186  CALL check_allocate(ierr,' IKLBOR')
1187  ALLOCATE(ifabor(nbtet,4),stat=ierr)
1188  CALL check_allocate(ierr,' IFABOR')
1189 
1190 ! BEGIN MODIF VSTOBIAC
1191  IF (.NOT. format_med) THEN
1192 
1193  ALLOCATE(ikle(nbtet,4),stat=ierr)
1194  CALL check_allocate(ierr,' IKLE')
1195 
1196  DO ielem = 1, nbtet
1197  DO i = 1,4
1198  ikle(ielem,i ) = ikles(ielem, i)
1199  END DO
1200  END DO
1201 
1202  DEALLOCATE(ikles)
1203 
1204  ENDIF
1205 !
1206  IF (nelin .GT. 0) THEN
1207  ALLOCATE(iklein(nelin,4),stat=ierr)
1208  CALL check_allocate(ierr,' IKLEIN')
1209 !
1210  DO ielem = 1, nelin
1211  DO i = 1,4
1212  iklein(ielem,i ) = iklestrin(ielem, i)
1213  END DO
1214  END DO
1215  ELSE
1216  ALLOCATE(iklein(1,4),stat=ierr)
1217  ENDIF
1218  DEALLOCATE(iklestrin)
1219 !
1220  WRITE(lu,*) 'FIN DE LA COPIE DE LA CONNECTIVITE INITIALE'
1221 !
1222  ALLOCATE(nbor(nptfr),stat=ierr)
1223  CALL check_allocate(ierr,' NBOR')
1224 !
1225  DO ielem = 1, nptfr
1226  nbor(ielem) = nbor2(ielem)
1227  ENDDO
1228 !
1229  DEALLOCATE(nbor2)
1230 !
1231  WRITE(lu,*) 'PARTEL_VOISIN31'
1232 !
1233  CALL voisin31(ifabor,nbtet,nbtet,31,
1234  & ikle,nbtet,npoint,nbor,nptfr,
1235  & lihbor,2,indpu,iklestri,nbtri)
1236 !
1237  WRITE(lu,*) 'FIN DE PARTEL_VOISIN31'
1238 !
1239  ALLOCATE(lihbor(nptfr),stat=ierr)
1240  CALL check_allocate(ierr,'LIHBOR')
1241  ALLOCATE(nulone(nbtri,3),stat=ierr)
1242  CALL check_allocate(ierr,' NULONE')
1243 !
1244  CALL elebd31( nelbor, nulone, iklbor,
1245  & ifabor, nbor, ikle,
1246  & nbtet, nbtri, nbtet, npoint,
1247  & nptfr,31)
1248 !
1249  DEALLOCATE(lihbor)
1250  DEALLOCATE(nulone)
1251 !
1252  WRITE(lu,*) 'FIN DE PARTEL_ELEBD31'
1253  ALLOCATE(number_tria(npoint),stat=ierr)
1254  CALL check_allocate(ierr,'NUMBER_TRIA')
1255  number_tria = 0
1256 !
1257  max_tria=0
1258  DO j = 1, nbtri
1259  k = 3*(j-1)+1
1260  ikle1 = iklestri(k)
1261  ikle2 = iklestri(k+1)
1262  ikle3 = iklestri(k+2)
1263  the_tri=ikle1
1264  IF (ikle2 < the_tri) the_tri=ikle2
1265  IF (ikle3< the_tri) the_tri=ikle3
1266  number_tria(the_tri)=number_tria(the_tri)+1
1267  END DO
1268  max_tria=maxval(number_tria)
1269 !
1270  DEALLOCATE(number_tria)
1271 !
1272  ALLOCATE(tri_ref(npoint,0:max_tria),stat=ierr)
1273  CALL check_allocate(ierr,' TRI_REF')
1274  tri_ref=0
1275  DO j = 1, nbtri
1276  k = 3*(j-1)+1
1277  ikle1 = iklestri(k)
1278  ikle2 = iklestri(k+1)
1279  ikle3 = iklestri(k+2)
1280  the_tri=ikle1
1281  IF (ikle2 < the_tri) the_tri=ikle2
1282  IF (ikle3< the_tri) the_tri=ikle3
1283  tri_ref(the_tri,0)=tri_ref(the_tri,0)+1
1284  pos=tri_ref(the_tri,0)
1285  tri_ref(the_tri,pos)=j
1286  END DO
1287 
1288  ALLOCATE(tettri(4*nbtet),stat=ierr)
1289  CALL check_allocate(ierr,' TETTRI')
1290  ALLOCATE(tettri2(nbtet),stat=ierr)
1291  CALL check_allocate(ierr,' TETTRI2')
1292  tettri(:) =-1
1293  tettri2(:) =0
1294 
1295  DO ieleb = 1,nbtri
1296  ielem = nelbor(ieleb)
1297  ikle1 = nbor(iklbor(ieleb,1))
1298  ikle2 = nbor(iklbor(ieleb,2))
1299  ikle3 = nbor(iklbor(ieleb,3))
1300  the_tri=ikle1
1301  IF (ikle2 < the_tri) the_tri=ikle2
1302  IF (ikle3<the_tri) the_tri=ikle3
1303  pos=tri_ref(the_tri,0)
1304  is = .false.
1305  m = -1
1306  DO jj = 1, pos
1307  j=tri_ref(the_tri,jj)
1308  k = 3*(j-1)+1
1309  IF ((ikle1.EQ.iklestri(k)).AND.
1310  & (ikle2.EQ.iklestri(k+1)).AND.
1311  & (ikle3.EQ.iklestri(k+2))) THEN
1312  is = .true.
1313  ELSE IF ((ikle1.EQ.iklestri(k)).AND.
1314  & (ikle3.EQ.iklestri(k+1)).AND.
1315  & (ikle2.EQ.iklestri(k+2))) THEN
1316  is = .true.
1317  ELSE IF ((ikle2.EQ.iklestri(k)).AND.
1318  & (ikle1.EQ.iklestri(k+1)).AND.
1319  & (ikle3.EQ.iklestri(k+2))) THEN
1320  is = .true.
1321  ELSE IF ((ikle2.EQ.iklestri(k)).AND.
1322  & (ikle3.EQ.iklestri(k+1)).AND.
1323  & (ikle1.EQ.iklestri(k+2))) THEN
1324  is = .true.
1325  ELSE IF ((ikle3.EQ.iklestri(k)).AND.
1326  & (ikle1.EQ.iklestri(k+1)).AND.
1327  & (ikle2.EQ.iklestri(k+2))) THEN
1328  is = .true.
1329  ELSE IF ((ikle3.EQ.iklestri(k)).AND.
1330  & (ikle2.EQ.iklestri(k+1)).AND.
1331  & (ikle1.EQ.iklestri(k+2))) THEN
1332  is = .true.
1333  ENDIF
1334  IF (is) THEN
1335  m = j
1336  EXIT
1337  ENDIF
1338  ENDDO
1339  DO i = 1,4
1340  IF (ifabor(ielem,i).LE.0) THEN
1341  IF ((ikle1.EQ.(ikle(nelbor(ieleb),somfac(1,i))))
1342  & .AND.(ikle2.EQ.(ikle(nelbor(ieleb),somfac(2,i))))
1343  & .AND. (ikle3.EQ.(ikle(nelbor(ieleb),somfac(3,i)))))
1344  & THEN
1345  ni = tettri2(ielem)
1346  n = 4*(ielem-1)+ni+1
1347  tettri(n) = m
1348 ! WRITE(*,*) N, '---> ',M
1349  tettri2(ielem) = ni + 1
1350  ENDIF
1351  ENDIF
1352  END DO
1353  ENDDO
1354 !
1355  DEALLOCATE(ifabor)
1356  DEALLOCATE(iklbor)
1357  DEALLOCATE(nelbor)
1358  DEALLOCATE(tri_ref)
1359  DEALLOCATE(ikle)
1360 !
1361  CALL system_clock(count=temps_sc(3),count_rate=parsec)
1362 !------------------------------------------------------- FIN VERSION F.D
1363 
1364 ! 3B. CONSTRUCTION DE NODES1/NODES2/NODES3: CONNECTIVITE INVERSE NOEUD > TETRA
1365 ! POUR L'ECRITURE A LA VOLEE DES UNV LOCAUX
1366 !---------------
1367 ! PARCOURS DES MAILLES POUR CONNAITRE LE NOMBRE DE MAILLES QUI
1368 ! LES REFERENCE
1369  ALLOCATE(nodes1(npoint),stat=ierr)
1370  CALL check_allocate(ierr,' NODES1')
1371  nodes1(:)=0
1372  DO i=1,nbtet
1373  DO k=1,4
1374  ikleb=iklestet(4*(i-1)+k)
1375  nodes1(ikleb)=nodes1(ikleb)+1
1376  ENDDO
1377  ENDDO
1378 ! NOMBRE DE REFERENCEMENT DE POINTS ET POINTEUR NODES2 VERS NODES3
1379 ! LE IEME POINT A SA LISTE DE TETRA (EN NUMEROTATION LOCALE TETRA)
1380 ! DE NODES3(NODES2(I)) A NODES3(NODES2(I)+NODES1(I)-1)
1381  ALLOCATE(nodes2(npoint+1),stat=ierr)
1382  CALL check_allocate(ierr,' NODES2')
1383  compt=0
1384  nodes2(1)=1
1385  DO i=1,npoint
1386  compt=compt+nodes1(i)
1387  nodes2(i+1)=compt+1
1388  ENDDO
1389 ! POUR UN NOEUDS DONNE, QU'ELLES SONT LES MAILLES QUI LE CONCERNENT
1390  ALLOCATE(nodes3(compt),stat=ierr)
1391  CALL check_allocate(ierr,' NODES3')
1392  nodes3(:)=-1
1393  DO i=1,nbtet
1394  DO k=1,4
1395  ikleb=iklestet(4*(i-1)+k)
1396  ni=nodes2(ikleb)
1397  nf=ni+nodes1(ikleb)-1
1398  nt=-999
1399  DO n=ni,nf ! ON CHERCHE LE PREMIER INDICE DE LIBRE DE NODES3
1400  IF (nodes3(n)==-1) THEN
1401  nt=n
1402  EXIT
1403  ENDIF
1404  ENDDO ! EN N
1405  IF (nt==-999) THEN
1406  GOTO 146 ! PB DE DIMENSIONNEMENT DE VECTEURS NODESI
1407  ELSE
1408  nodes3(nt)=i ! NUMERO LOCAL DU TETRA I ASSOCIE AU NOEUD NT
1409  ENDIF
1410  ENDDO
1411  ENDDO
1412 
1413 ! 3C. CONSTRUCTION DE NODES1T/NODES2T/NODES3T: CONNECTIVITE INVERSE NOEUD > TRIA
1414 ! POUR LA COULEUR DES NOEUDS (DIRICHLET SUR L'INTERFACE)
1415 !---------------
1416  ALLOCATE(nodes1t(npoint),stat=ierr)
1417  CALL check_allocate(ierr,' NODES1T')
1418  nodes1t(:)=0
1419  DO i=1,nbtri
1420  DO k=1,3
1421  ikleb=iklestri(3*(i-1)+k)
1422  nodes1t(ikleb)=nodes1t(ikleb)+1
1423  ENDDO
1424  ENDDO
1425 
1426  ALLOCATE(nodes2t(npoint+1),stat=ierr)
1427  CALL check_allocate(ierr,' NODES2T')
1428  compt=0
1429  nodes2t(1)=1
1430  DO i=1,npoint
1431  compt=compt+nodes1t(i)
1432  nodes2t(i+1)=compt+1
1433  ENDDO
1434  ALLOCATE(nodes3t(compt),stat=ierr)
1435  CALL check_allocate(ierr,' NODES3T')
1436  nodes3t(:)=-1
1437  DO i=1,nbtri
1438  DO k=1,3
1439  ikleb=iklestri(3*(i-1)+k)
1440  ni=nodes2t(ikleb)
1441  nf=ni+nodes1t(ikleb)-1
1442  nt=-999
1443  DO n=ni,nf ! ON CHERCHE LE PREMIER INDICE DE LIBRE DE NODES3T
1444  IF (nodes3t(n)==-1) THEN
1445  nt=n
1446  EXIT
1447  ENDIF
1448  ENDDO ! EN N
1449  IF (nt==-999) THEN
1450  GOTO 146 ! PB DE DIMENSIONNEMENT DE VECTEURS NODESI
1451  ELSE
1452  nodes3t(nt)=i ! NUMERO LOCAL DU TETRA I ASSOCIE AU NOEUD NT
1453  ENDIF
1454  ENDDO
1455  ENDDO
1456  CALL system_clock(count=temps_sc(4),count_rate=parsec)
1457  WRITE(lu,*)' TEMPS CONNECTIVITE INVERSE PART1/ PART2',
1458  & (1.0*(temps_sc(3)-temps_sc(2)))/(1.0*parsec),'/',
1459  & (1.0*(temps_sc(4)-temps_sc(3)))/(1.0*parsec),' SECONDS'
1460 
1461 !----------------------------------------------------------------------
1462 ! 4. PARTITIONING
1463 !---------------
1464 
1465 ! DO I=1,4*NBTET
1466 ! WRITE(LU,*) 'TETTRIALPHA',TETTRI(I)
1467 ! ENDDO
1468 ! DO I=1,NBTET
1469 ! WRITE(LU,*) 'TETTRIBETA',TETTRI2(I)
1470 ! ENDDO
1471 !
1472 !----------------------------------------------------------------------
1473 ! NEW METIS INTERFACE (>= VERSION 5) :
1474 !
1475  ALLOCATE(epart(nbtet),stat=ierr)
1476  CALL check_allocate (ierr, 'EPART')
1477  ALLOCATE (npart(npoint),stat=ierr)
1478  CALL check_allocate (ierr, 'NPART')
1479 !
1480 !----------------------------------------------------------------------
1481 ! CALL METIS : MESH PARTITIONNING
1482 !
1483  CALL system_clock(count=temps_sc(5),count_rate=parsec)
1484 !
1485  WRITE(lu,*)' '
1486  WRITE(lu,*)' STARTING METIS MESH PARTITIONING------------------+'
1487 !
1488  CALL partitioner(pmethod, nbtet, npoint, 4, nparts, iklestet,
1489  & epart, npart)
1490 !
1491  CALL system_clock(count=temps_sc(6),count_rate=parsec)
1492 !
1493  WRITE(lu,*)' '
1494  WRITE(lu,*)' END METIS MESH PARTITIONING------------------+'
1495  WRITE(lu,*)' TEMPS CONSOMME PAR METIS ',
1496  & (1.0*(temps_sc(6)-temps_sc(5)))/(1.0*parsec),' SECONDS'
1497  WRITE(lu,80) nelemtotal,npoint
1498  WRITE(lu,81) nbtet,nbtri
1499  WRITE(lu,82) nparts
1500  WRITE(lu,*) 'SORTIE DE METIS CORRECTE'
1501 !
1502 !
1503 !D ******************************************************
1504 !D LOOP OVER THE TETRA TO COMPUTER THE NUMBER AND THE LABEL
1505 !D OF FINITE ELEMENTS ASSIGNED TO EACH SUBDOMAIN
1506 !D ******************************************************
1507 !D COMPUTATION OF THE MAXIMUM NUMBER OF FINITE ELEMENTS ASSIGNED TO ONE SUBDOMAIN
1508  nelem_p(:)=0
1509  npoin_p(:)=0
1510  DO i=1,nbtet
1511  nelem_p(epart(i))=nelem_p(epart(i))+1
1512  END DO
1513  max_nelem_p=maxval(nelem_p)
1514  WRITE(lu,*) 'NB MAX OF TETRAS PER SUBDOMAIN : ',max_nelem_p
1515  WRITE(lu,*) 'NB OF TETRA PER SUBDOMAIN :'
1516  DO i=1,nparts
1517  WRITE(lu,*) i, nelem_p(i)
1518  ENDDO
1519  nelem_p(:)=0
1520 !D ALLOCATION OF THE ELEGL ARRAY
1521  ALLOCATE(elegl(max_nelem_p,nparts),stat=ierr)
1522 !D ELEGL IS THE FILLED
1523  CALL check_allocate(ierr,'ELEGL')
1524  DO i=1,nbtet
1525  nelem_p(epart(i))=nelem_p(epart(i))+1
1526  elegl(nelem_p(epart(i)),epart(i))=i
1527  END DO
1528 !D COMPUTE THE MAXIMUM OF NODES ASSIGNED TO ONE SUBDOMAIN
1529 
1530  ALLOCATE(nodelg(npoint,nparts),stat=ierr)
1531  CALL check_allocate(ierr,'NODELG')
1532 
1533  nodelg(:,:)=0
1534 !D FOR EACH SUBDOMAIN IDD
1535  DO idd=1,nparts
1536 !D LOOP ON THE FINITE ELEMENTS IELEM ASSIGNED TO SUBDOMAIN IDD
1537  DO pos=1,nelem_p(idd)
1538  ielem=elegl(pos,idd)
1539  n=4*(ielem-1)+1
1540 !D LOOP OF THE NODE CONTAINED IN IELEM
1541  DO k=0,3
1542  node=iklestet(n+k)
1543  IF (nodelg(node,idd) .EQ. 0) THEN
1544  npoin_p(idd)=npoin_p(idd)+1
1545  nodelg(node,idd)=npoin_p(idd)
1546  END IF
1547  END DO
1548  END DO
1549  END DO
1550 !D ALLOCATION AND FILLING OF THE NODEGL ARRAY
1551  max_npoin_p=maxval(npoin_p)
1552 
1553  WRITE(lu,*) 'NB MAX OF POINT PER SUBDOMAIN :', max_npoin_p
1554  WRITE(lu,*) 'NB OF POINT PER SUBDOMAIN :'
1555  DO i=1,nparts
1556  WRITE(lu,*) i, npoin_p(i)
1557  ENDDO
1558 !
1559  ALLOCATE(nodegl(max_npoin_p,nparts),stat=ierr)
1560  CALL check_allocate(ierr,'NODEGL')
1561  nodegl(:,:)=0
1562  DO idd=1,nparts
1563  DO node=1,npoint
1564  IF (nodelg(node,idd) .NE. 0) THEN
1565  nodegl(nodelg(node,idd),idd)=node
1566  END IF
1567  END DO
1568  END DO
1569 !
1570 !----------------------------------------------------------------------
1571 ! 5A. ALLOCATIONS POUR ECRITURE DES FICHIERS .UNV/.LOG ASSOCIANT UN SOUS-DOMAINE
1572 ! PAR PROC
1573 !------------
1574 
1575  nameinp2=nameinp
1576  namelog2=namelog
1577  blanc=' '
1578  moins1='-1'
1579  ALLOCATE(nodes4(npoint),stat=ierr)
1580  CALL check_allocate(ierr,' NODES4')
1581 !$$$ NODES4(:)=-1
1582  ALLOCATE(knolg(npoint),stat=ierr) ! C'EST SOUS-OPTIMAL EN
1583  CALL check_allocate(ierr,' KNOLG')! TERME DE DIMENSIONNEMENT
1584  knolg(:)=-1 ! MAIS PLUS RAPIDE POUR LE REMPLISSAGE ULTERIEUR
1585 !
1586 ! PARAMETRE NBSDOMVOIS (NOMBRE DE SOUS DOMAINES VOISINS+2)
1587 !
1588  ALLOCATE(nachb(nbsdomvois,npoint),stat=ierr)
1589  CALL check_allocate(ierr,' NACHB')
1590  nachb(1,:)=0
1591  DO j=2,nbsdomvois-1
1592  nachb(j,:)=-1
1593  ENDDO
1594  ALLOCATE(triunv(4*nbtri),stat=ierr)
1595  CALL check_allocate(ierr, 'TRIUNV')
1596 !
1597 !
1598 ! 5B. RECHERCHE DE LA VRAI COULEUR AUX NOEUDS POUR EVITER LES PBS DE DIRICHLET
1599 ! AUX INTERFACES
1600 !---------------
1601  ncolor2(:)=-1
1602  DO j=1,npoint ! BOUCLE SUR TOUS LES POINTS DU MAILLAGES
1603  ni=nodes2t(j)
1604  nf=ni+nodes1t(j)-1
1605 
1606  DO n=ni,nf ! BOUCLE SUR LES TETRA CONTENANT LE POINT J
1607  numtet=nodes3t(n) ! TRIA DE NUMERO LOCAL NUMTET
1608  numtrig=convtri(numtet) ! NUMERO GLOBAL DU TRIANGLE
1609  color1=ecolor(numtrig) ! COULEUR DU NOEUD AVEC CE TRIA
1610  color2=ncolor2(j)
1611 
1612  IF (color2 > 0) THEN ! ON PRIORISE LES COULEURS
1613  pr1=0
1614  pr2=0
1615  DO l=1,nbcolor
1616  IF (priority(l)==color1) THEN
1617  pr1=l
1618  ENDIF
1619  IF (priority(l)==color2) THEN
1620  pr2=l
1621  ENDIF
1622  ENDDO
1623  IF ((pr1==0).OR.(pr2==0)) GOTO 154
1624  IF (pr1<pr2) ncolor2(j)=color1 ! ON CHANGE DE COULEUR
1625  ELSE ! PREMIERE FOIS QUE CE NOEUD EST TRAITE
1626  ncolor2(j)=color1
1627  ENDIF
1628  ENDDO
1629  ENDDO
1630 
1631  CALL system_clock(count=temps_sc(7),count_rate=parsec)
1632 
1633 ! DO IELEM = 1, NPOINT
1634 ! WRITE(LU,*) 'NCOLOR2',NCOLOR2(IELEM)
1635 ! ENDDO
1636 
1637 ! DO IELEM = 1, NBCOLOR
1638 ! WRITE(LU,*) 'PRIOR',PRIORITY(IELEM)
1639 ! ENDDO
1640 !
1641 ! OB D
1642 !--------------
1643 ! RAJOUT POUR TENIR COMPTE DES COULEURS DES NOEUDS DE TETRAS LIES
1644 ! AU TRIA DE BORD ET SITUES DANS D'AUTRES SD
1645 !--------------
1646 !
1647  ALLOCATE(tetcolor(nbtet,4),stat=ierr)
1648  CALL check_allocate(ierr,' TETCOLOR')
1649  tetcolor(:,:)=.false.
1650  nbretouche=0
1651  DO iptfr=1,nptfr ! BOUCLE SUR TOUS LES POINTS DE BORD
1652  j=nbor(iptfr)
1653 ! ON NE FAIT QQE CHOSE (EVENTUELLEMENT) QUE SI IL Y A UN TRIA
1654 ! DE BORD (ECOLOR>0 ET NCOLOR2 !=-1). GRACE AU TRAITEMENT PRECEDENT
1655 ! ON S'EN REND COMPTE DIRECTEMENT VIA NCOLOR2.
1656  linter=.false.
1657  nbtetj=nodes1(j) ! NBRE DE TETRA RATTACHES A CE NOEUD
1658  ni=nodes2(j) ! ADRESSE DANS NODES3 DU PREMIER
1659  nf=ni+nbtetj-1
1660  IF (ncolor2(j) > 0) THEN
1661 ! ON CHERCHE A SAVOIR SI LE NOEUD EST A L'INTERFACE LINTER=.TRUE.
1662  DO n=ni,nf ! BOUCLE SUR LES TETRA CONTENANT LE POINT J
1663  nt=nodes3(n) ! TETRA DE NUMERO LOCAL NT
1664  IF (n.EQ.ni) THEN
1665  iddnt=epart(nt)
1666  ELSE
1667  IF (epart(nt) /= iddnt) THEN
1668  linter=.true.
1669  GOTO 20 ! ON A LE RENSEIGNEMENT DEMANDE, ON SORT
1670  ENDIF
1671  ENDIF
1672  ENDDO ! FIN BOUCLE SUR LES TETRAS
1673  20 CONTINUE
1674 ! LE NOEUD J EST UN NOEUD D'INTERFACE. ON VA COMMUNIQUER AU NOEUD
1675 ! CORRESPONDANT DES TETRAS (SI UN TRIA DE BORD N'EST PAS SUR CETTE
1676 ! FACE AUXQUEL CAS LE PB EST DEJA REGLE), LA BONNE COULEUR.
1677  IF (linter) THEN
1678  DO n=ni,nf ! BOUCLE SUR LES TETRA CONTENANT LE POINT J
1679  nt=nodes3(n) ! TETRA DE NUMERO LOCAL NT
1680 ! ON VA TRIER LES CAS NON PATHOLOGIQUES ET TRES COURANT DE TETRA
1681 ! DONT UNE FACE COINCIDE AVEC CE TRIANGLE
1682  IF (tettri2(nt)>0) THEN !TETRA CONCERNE PAR UN TRIA
1683  nit=4*(nt-1)+1
1684  nft=nit+tettri2(nt)-1
1685  DO mt=nit,nft ! BOUCLE SUR LES TRIA DU TETRA
1686  numtri=tettri(mt) ! NUM LOCAL DU TRIA
1687  numtrib=3*(numtri-1)+1
1688  ikle1=iklestri(numtrib) ! NUMERO GLOBAUX DES NOEUDS DU TRIA
1689  ikle2=iklestri(numtrib+1)
1690  ikle3=iklestri(numtrib+2)
1691 ! CE POINT J APPARTIENT DEJA A UN TRIA ACOLLE AU TETRA
1692 ! ON SAUTE LE TETRA NT
1693  IF ((ikle1==j).OR.(ikle2==j).OR.(ikle3==j)) THEN
1694 ! POUR TESTS
1695 ! WRITE(LU,*)'JE SAUTE LE TETRA ',NT,EPART(NT),
1696 ! & TETTRI2(NT),' NODES ',J
1697  GOTO 21
1698  ENDIF
1699  ENDDO
1700  ENDIF ! FIN SI TETTRI
1701 ! LE TETRA NT EST POTENTIELLEMENT OUBLIE, ON LE TRAITE AU CAS OU
1702 ! LE PARTAGE SE FERA DANS ESTEL3D/READ_CONNECTIVITY
1703  numtetb=4*(nt-1)+1
1704  DO l=1,4
1705  ikle1=iklestet(numtetb+l-1) ! NUMERO GLOBAUX DES NOEUDS DU TETRA
1706  IF (ikle1==j) THEN
1707  tetcolor(nt,l)=(tetcolor(nt,l).OR..true.)
1708  nbretouche=nbretouche+1
1709  ENDIF
1710  ENDDO ! EN L
1711  21 CONTINUE
1712  ENDDO ! FIN BOUCLE SUR LES TETRAS
1713  ENDIF ! FIN SI LINTER
1714  ENDIF ! FIN SI SUR NCOLOR2
1715  ENDDO ! FIN BOUCLE SUR LES POINTS DE BORD
1716 ! OB F
1717  CALL system_clock(count=temps_sc(8),count_rate=parsec)
1718  WRITE(lu,*)' NOMBRE DE RETOUCHE DU PARTITIONNEMENT (PART2): ',
1719  & nbretouche
1720  WRITE(lu,*)' TEMPS DE RETOUCHE DU PARTITIONNEMENT PART1/PART2',
1721  & (1.0*(temps_sc(7)-temps_sc(6)))/(1.0*parsec),'/',
1722  & (1.0*(temps_sc(8)-temps_sc(7)))/(1.0*parsec),' SECONDS'
1723 !$$$ WRITE(LU,*)'IDEM VERSION DE REFERENCE'
1724 
1725 ! 5C. REMPLISSAGE EFFECTIF DU UNV PAR SD
1726 !---------------
1727  ibid = 1
1728 !
1729  IF (nelin .GT. 0) THEN
1730  ALLOCATE(deja_trouve(nelin),stat=ierr)
1731  ELSE
1732  ALLOCATE(deja_trouve(1),stat=ierr)
1733  ENDIF
1734  CALL check_allocate(ierr,'DEJA_TROUVE')
1735  deja_trouve(:)=.false.
1736 !
1737 ! INITIALISATION OF THE SIZE OF THE FILE'S NAME
1738  i_leninp = len_trim(nameinp2)
1739  i_lenlog = len_trim(namelog2)
1740 !
1741  DO idd=1,nparts ! BOUCLE SUR LES SOUS-DOMAINES
1742 
1743 ! NOMBRE DE TRIANGLES POUR CE SOUS-DOMAINE
1744  nbtriidd=0
1745 ! NOM DU FICHIER UNV PAR SOUS-DOMAINE
1746  nameinp2(i_leninp+1:i_leninp+11) = extens(nparts-1,idd-1)
1747  OPEN(ninp2,file=nameinp2,status='UNKNOWN',form='FORMATTED',
1748  & err=132)
1749  rewind(ninp2)
1750 
1751 ! NOM DU FICHIER LOG PAR SOUS-DOMAINE
1752  namelog2(i_lenlog+1:i_lenlog+11) = extens(nparts-1,idd-1)
1753  OPEN(nlog2,file=namelog2,status='UNKNOWN',form='FORMATTED',
1754  & err=133)
1755  rewind(nlog2)
1756 
1757 ! TITRE (UNV PAR SD)
1758  WRITE(ninp2,60,err=112)blanc,moins1
1759  WRITE(ninp2,61,err=112)nsec1
1760  WRITE(ninp2,62,err=112)titre
1761  titre = ' '
1762  WRITE(ninp2,62,err=112)titre
1763  WRITE(ninp2,62,err=112)titre
1764  WRITE(ninp2,62,err=112)titre
1765  WRITE(ninp2,62,err=112)titre
1766  WRITE(ninp2,62,err=112)titre
1767  WRITE(ninp2,62,err=112)titre
1768  WRITE(ninp2,60,err=112)blanc,moins1
1769 !
1770 ! BLOC SUR LES COORDONNEES/COULEURS DES NOEUDS (UNV PAR SD)
1771  WRITE(ninp2,60,err=112)blanc,moins1
1772  WRITE(ninp2,61,err=112)nsec2
1773  compt=1
1774  nodes4(:)=-1
1775 !D NEW VERSION OF THE LOOP TO REDUCE THE COMPUTING TIME
1776  DO pos_node=1,npoin_p(idd) ! BOUCLE SUR TOUS LES POINTS DU MAILLAGES
1777  j=nodegl(pos_node,idd)
1778 !D PREVIOUS VERSION OF THE LOOP
1779 !D NI=NODES2(J)
1780 !D NF=NI+NODES1(J)-1
1781 !D DO N=NI,NF ! BOUCLE SUR LES TETRA CONTENANT LE POINT J
1782 !D NT=NODES3(N) ! TETRA DE NUMERO LOCAL NT
1783 !D IF (EPART(NT)==IDD) THEN ! C'EST UNE MAILLE DU SOUS-DOMAINE
1784  WRITE(ninp2,63,err=112)compt,ibid,ibid,ncolor2(j)
1785  WRITE(ninp2,64,err=112)x1(j),y1(j),z1(j)
1786  nodes4(j)=compt ! LE NOEUD J A LE NUMERO COMPT
1787  ! POUR LE SOUS-DOMAINE IDD
1788 ! POUR PARALLELISME TELEMAC
1789  knolg(compt)=j ! CONVERSION SD (LOCAL)-->MAILLAGE ENTIER (GLOBAL)
1790  k=nachb(1,j) ! NBRE DE SD CONTENANT LE NOEUD J
1791  nachblog=.true.
1792  DO l=1,k ! NOEUD DEJA CONCERNE PAR CE SD ?
1793  IF (nachb(1+l,j)==idd) nachblog=.false. ! OUI
1794  ENDDO
1795  IF (nachblog) THEN ! NON
1796  k=nachb(1,j)+1
1797  IF (k.GT.nbsdomvois-2) GOTO 151
1798  nachb(k+1,j)=idd ! NOEUD GLOBAL J CONCERNE PAR IDD
1799  nachb(1,j)=k ! SA MULTIPLICITE
1800  ENDIF
1801  compt=compt+1
1802 ! GOTO 10 ! ON PASSE AU NOEUD SUIVANT
1803 ! ENDIF ! EN EPART
1804 ! ENDDO ! EN N
1805 ! 10 CONTINUE
1806  ENDDO ! EN J
1807 ! POUR TESTS
1808 ! DO I=1,NPOINT
1809 ! WRITE(LU,*)'GLOBAL NUMERO POINT: ',I,' LOCAL: ',NODES4(I)
1810 ! ENDDO
1811  npointsd(idd)=compt-1 ! NOMBRE DE NOEUDS DU SOUS-DOMAINE IDD
1812  WRITE(ninp2,60,err=112)blanc,moins1
1813 
1814 ! BLOC SUR LES CONNECTIVITES/COULEURS DES MAILLES (UNV PAR SD)
1815  WRITE(ninp2,60,err=112)blanc,moins1
1816  WRITE(ninp2,61,err=112)nsec3
1817  compt=1
1818  ibid = 1
1819 !D PREVIOUS VERSION OF THE LOOP
1820 !D DO J=1,NELEMTOTAL
1821 !D IF (TYPELEM(J,1)==111) THEN ! C'EST UN TETRAEDRE
1822 !D NUMTET=TYPELEM(J,2) ! NUM LOCAL DU TETRA DANS LA LISTE DES TETRAS
1823 !D IF (EPART(NUMTET)==IDD) THEN
1824 
1825  DO pos=1,nelem_p(idd)
1826  ! BOUCLE SUR TETRA ET TRIA POUR ECOLOR
1827  j=elegl(pos,idd)
1828  numtet=typelem(j,2) ! NUM LOCAL DU TETRA DANS LA LISTE DES TETRAS
1829  elem=111
1830 ! OB D
1831 ! PRETRAITEMENT POUR LES EVENTUELS PB DE COULEURS DES NOEUDS DE TETRAS
1832 ! A L'INTERFACE
1833  ibidc=0
1834  IF (tetcolor(numtet,1)) ibidc=ibidc+1000
1835  IF (tetcolor(numtet,2)) ibidc=ibidc+ 200
1836  IF (tetcolor(numtet,3)) ibidc=ibidc+ 30
1837  IF (tetcolor(numtet,4)) ibidc=ibidc+ 4
1838 ! POUR MONITORING
1839 ! IF (IBIDC/=0) WRITE(LU,*)'IDD',IDD,'PARTEL',J,COMPT,IBIDC
1840 ! IDEM VERSION DE REFERENCE
1841 ! IBIDC=0
1842 ! OB F
1843  WRITE(ninp2,65,err=112)compt,elem,-ibidc,ibid,ecolor(j),4
1844  IF (ecolor(j).LE.0) print*,'PB WRITE COLOR',j,ecolor(j)
1845  compt=compt+1
1846  n=4*(numtet-1)+1
1847  ikle1=nodes4(iklestet(n))
1848  ikle2=nodes4(iklestet(n+1))
1849  ikle3=nodes4(iklestet(n+2))
1850  ikle4=nodes4(iklestet(n+3))
1851  WRITE(ninp2,66,err=112)ikle1,ikle2,ikle3,ikle4
1852  IF ((ikle1.LT.0).OR.(ikle2.LT.0).OR.(ikle3.LT.0)
1853  & .OR.(ikle4.LT.0)) GOTO 147
1854  IF (tettri2(numtet).NE.0) THEN
1855  ni=4*(numtet-1)+1
1856  nf=ni+tettri2(numtet)-1
1857  DO m=ni,nf ! ON PARCOURT LES TRIANGLES DE BORD ASSOCIES
1858  numtri=tettri(m) ! AU NUMTET TETRAEDRE; NUM LOCAL DU TRIA
1859  numtrig=convtri(numtri) ! NUMERO GLOBAL DU TRIANGLE
1860  elem=91
1861  triunv(4*nbtriidd+1)=ecolor(numtrig)
1862  n=3*(numtri-1)+1
1863  ikle1=nodes4(iklestri(n))
1864  ikle2=nodes4(iklestri(n+1))
1865  ikle3=nodes4(iklestri(n+2))
1866  triunv(4*nbtriidd+2)=ikle1
1867  triunv(4*nbtriidd+3)=ikle2
1868  triunv(4*nbtriidd+4)=ikle3
1869  nbtriidd=nbtriidd+1
1870 !
1871  IF ((ikle1.LT.0).OR.(ikle2.LT.0).OR.(ikle3.LT.0))
1872  & GOTO 147
1873 !
1874  ENDDO ! EN M
1875  ENDIF ! EN TETTRI2
1876  ENDDO ! EN J
1877 
1878 ! MAINTENANT ON PEUX RECOPIER LE BLOC DES TRIANGLES !
1879  elem=91
1880  DO j=1,nbtriidd
1881  WRITE(ninp2,65,err=112)compt,elem,ibid,ibid,
1882  & triunv(4*(j-1)+1),3
1883  ikle1=triunv(4*(j-1)+2)
1884  ikle2=triunv(4*(j-1)+3)
1885  ikle3=triunv(4*(j-1)+4)
1886  WRITE(ninp2,67,err=112)ikle1,ikle2,ikle3
1887  compt=compt+1
1888  ENDDO ! EN J
1889 !
1890  elem=91
1891 ! BOUCLE SURDIMENSIONNEE, ON BOUCLE SUR LE NOMBRE DE SURFACE INTERNE DU MAILLAGE GLOBAL...
1892  IF (nelin .GT. 0) THEN
1893  DO j=1,nelin
1894  IF (deja_trouve(j)) cycle
1895  ikle1=nodes4(iklein(j,2))
1896  ikle2=nodes4(iklein(j,3))
1897  ikle3=nodes4(iklein(j,4))
1898  IF ((ikle1.EQ.-1).OR.(ikle2.EQ.-1).OR.(ikle3.EQ.-1)) cycle
1899  !
1900  ! MODIF STOBIAC
1901  found_tet = .false.
1902  !
1903  ! ON BOUCLE SUR LES TETRA DE CHAQUE NOEUD POUR VERIFIER QUE LES TROIS POINTS
1904  ! APPARTIENNENT A UN TETRAEDRE DE LA PARTITION
1905  !
1906  ! PASSAGE A LA NUMEROTATION GLOBALE
1907  ptri1 = nodegl(ikle1,idd)
1908  ptri2 = nodegl(ikle2,idd)
1909  ptri3 = nodegl(ikle3,idd)
1910  !
1911  deb1 = nodes2(ptri1)
1912  fin1 = deb1 + nodes1(ptri1)-1
1913  deb2 = nodes2(ptri2)
1914  fin2 = deb2 + nodes1(ptri2)-1
1915  deb3 = nodes2(ptri3)
1916  fin3 = deb3 + nodes1(ptri3)-1
1917  !
1918  DO ptet1 = deb1, fin1
1919  DO ptet2 = deb2, fin2
1920  IF (nodes3(ptet1).EQ.nodes3(ptet2)) THEN
1921  DO ptet3 = deb3, fin3
1922  IF (nodes3(ptet3).EQ.nodes3(ptet1)) THEN
1923  IF (epart(nodes3(ptet3)).EQ.idd) THEN
1924  found_tet = .true.
1925  ENDIF
1926  ENDIF
1927  IF (found_tet) EXIT
1928  ENDDO
1929  ENDIF
1930  IF (found_tet) EXIT
1931  ENDDO
1932  IF (found_tet) EXIT
1933  ENDDO
1934  !
1935  IF (.NOT.found_tet) cycle
1936  ! END MODIF STOBIAC
1937  !
1938  WRITE(ninp2,65,err=112) compt,elem,ibid,ibid,iklein(j,1),3
1939  WRITE(ninp2,67,err=112) ikle1,ikle2,ikle3
1940  compt = compt+1
1941  deja_trouve(j) = .true.
1942  ENDDO ! EN J
1943  ENDIF
1944 !
1945 !$$$ WRITE(LU,*) 'SUBDOMAIN',IDD,'INNERTRI',COMPT
1946 !
1947  WRITE(ninp2,60,err=112)blanc,moins1
1948 ! WRITE(NINP2,60,ERR=112)BLANC,MOINS1
1949 ! WRITE(NINP2,61,ERR=112)NSEC4
1950 ! WRITE(NINP2,68,ERR=112) 1,0,0,0,0,0,0,0
1951  CLOSE(ninp2)
1952  nelemsd(idd)=compt-1 ! NOMBRE DE MAILLES DU SOUS-DOMAINE IDD
1953 
1954 ! 5D. REMPLISSAGE EFFECTIF DU LOG PAR SD
1955 !---------------
1956 ! ELEMENT STANDARD DU FICHIER LOG (LOG PAR SD)
1957  WRITE(nlog2,51 ,err=113) npointsd(idd)
1958  WRITE(nlog2,52 ,err=113) nelemsd(idd)
1959  WRITE(nlog2,523,err=113) size_flux
1960 
1961 ! BEGIN MODIF V STOBIAC
1962 ! READ FAMILIES
1963 #if defined HAVE_MED
1964  IF (format_med) THEN
1965  WRITE(nlog2,53 ,err=113) nbfamily
1966  DO j=1,nbfamily+1
1967  WRITE(nlog2,50,err=113)'--'
1968  ENDDO
1969  ENDIF
1970 #endif
1971  IF (.NOT.format_med) THEN
1972  WRITE(nlog2,53 ,err=113) nbfamily-1
1973  DO j=1,nbfamily
1974  WRITE(nlog2,50,err=113)'--'
1975  ENDDO
1976  ENDIF
1977 ! END MODIF V STOBIAC
1978 
1979  ! ADDITION BY JP RENAUD ON 15/02/2007
1980  ! AS THE LIST OF PRIORITIES HAS MOVED IN ESTEL-3D FROM
1981  ! THE STEERING FILE TO THE LOG FILE, WE NEED TO WRITE "A"
1982  ! NUMBER OF EXTERNAL FACES + PRIORITY LIST HERE. AS THESE
1983  ! ARE NOT USED IN PARALLEL MODE, WE MERELY COPY THE LIST
1984  ! FROM THE ORIGINAL LOG FILE.
1985 
1986  WRITE(nlog2,531,err=113) nbcolor
1987  WRITE(unit=theformat,fmt=1000) nbcolor
1988 1000 FORMAT('(''PRIORITY :'',',i3,'(X,I3,))')
1989  theformat=trim(theformat)
1990 ! WRITE(LU,*) 'FORMATT =',THEFORMAT
1991  WRITE (nlog2,fmt=theformat(1:len(theformat)-1))
1992  & (priority(i), i=1, nbcolor)
1993 
1994  ! END ADDITION BY JP RENAUD
1995 
1996 ! KNOLG (LOG PAR SD)
1997  nt=npointsd(idd)
1998  ni=nt/6
1999  nf=nt-6*ni
2000  WRITE(nlog2,54,err=113)ni,nf
2001  DO j=1,ni
2002  WRITE(nlog2,540,err=113)(knolg(6*(j-1)+k),k=1,6)
2003  ENDDO
2004  IF (nf.EQ.1) THEN
2005  WRITE(nlog2,541,err=113)knolg(6*ni+1)
2006  ELSE IF (nf.EQ.2) THEN
2007  WRITE(nlog2,542,err=113)(knolg(6*ni+k),k=1,2)
2008  ELSE IF (nf.EQ.3) THEN
2009  WRITE(nlog2,543,err=113)(knolg(6*ni+k),k=1,3)
2010  ELSE IF (nf.EQ.4) THEN
2011  WRITE(nlog2,544,err=113)(knolg(6*ni+k),k=1,4)
2012  ELSE IF (nf.EQ.5) THEN
2013  WRITE(nlog2,545,err=113)(knolg(6*ni+k),k=1,5)
2014  ENDIF
2015  WRITE(nlog2,55,err=113)npoint ! NOMBRE DE NOEUD DU MAILLAGE
2016  ! INITIAL POUR ALLOCATION KNOGL DANS ESTEL
2017 !
2018  ENDDO ! BOUCLE SUR LES SOUS-DOMAINES
2019 
2020  DEALLOCATE(convtri)
2021  DEALLOCATE(typelem)
2022 
2023 ! 5E. TRAVAUX SUPPLEMENTAIRES POUR DETERMINER LE NACHB AVANT DE L'ECRIRE
2024 ! DANS LE LOG
2025 !---------------
2026  DO idd=1,nparts ! BOUCLE SUR LES SOUS-DOMAINES
2027 ! CONSTRUCTION ET DIMENSIONNEMENT DU NACHB PROPRE A CHAQUE SD
2028  compt=0
2029  nachb(nbsdomvois,:)=-1
2030  DO j=1,npoint ! BOUCLE SUR TOUS LES POINTS DU MAILLAGE
2031  n=nachb(1,j)
2032  IF (n>1) THEN ! POINT D'INTERFACE
2033  n=n+1
2034  DO k=2,n
2035  IF (nachb(k,j)==idd) THEN ! IL CONCERNE IDD
2036  compt=compt+1 ! "COMPT"IEME POINT D'INTERFACE DE IDD
2037  nachb(nbsdomvois,j)=compt ! A RETENIR COMME POINT D'INTERFACE
2038  ENDIF
2039  ENDDO ! FIN BOUCLE SUR LES SD DU POINT J
2040  ENDIF
2041  ENDDO ! FIN BOUCLE POINTS
2042  npointisd(idd)=compt ! NOMBRE DE POINTS D'INTERFACE DE IDD
2043 
2044 ! 5F. ON CONTINUE L'ECRITURE DU .LOG
2045 !-------------
2046  namelog2(i_lenlog+1:i_lenlog+11) = extens(nparts-1,idd-1)
2047  OPEN(nlog2,file=namelog2,status='OLD',form='FORMATTED',
2048  & position='APPEND',err=133)
2049  WRITE(nlog2,56,err=113) npointisd(idd)
2050  DO j=1,npoint
2051  IF (nachb(nbsdomvois,j)>0) THEN ! C'EST UN POINT D'INTERFACE DE IDD
2052  compt=0
2053  vectnb(:)=-1
2054  DO k=1,nbsdomvois-2 ! ON PREPARE L'INFO POUR LE NACHB TELEMAC
2055  IF (nachb(k+1,j)/= idd) THEN
2056  compt=compt+1
2057 ! ATTENTION A CELUI-CI, SUREMENT LIE AU NUMERO DE POINTS...
2058 ! OB D
2059  IF (compt.GT.nbsdomvois-3) GOTO 152
2060 ! OB F
2061  IF (nachb(k+1,j)>0) THEN
2062 ! ON STOCKE LE NUMERO DE PROC ET NON LE NUMERO DE SOUS-DOMAINE
2063 ! D'OU LA CONTRAINTE, UN PROC PAR SOUS-DOMAINE
2064  vectnb(compt)=nachb(k+1,j)-1
2065  ENDIF
2066  ENDIF
2067  ENDDO ! EN K
2068 ! WRITE(NLOG2,561,ERR=113)J,(VECTNB(K),K=1,NBSDOMVOIS-3)
2069  nt = nbsdomvois-3
2070  ni=nt/6
2071  nf=nt-6*ni+1
2072  WRITE(nlog2,640,err=113)nodelg(j,idd),(vectnb(k),k=1,5)
2073 ! WRITE(NLOG2,640,ERR=113)J,(VECTNB(K),K=1,5)
2074  DO l=2,ni
2075  WRITE(nlog2,640,err=113)(vectnb(6*(l-1)+k),k=0,5)
2076  ENDDO
2077  IF (nf.EQ.1) THEN
2078  WRITE(nlog2,641,err=113)vectnb(6*ni)
2079  ELSEIF (nf.EQ.2) THEN
2080  WRITE(nlog2,642,err=113)(vectnb(6*ni+k),k=0,1)
2081  ELSEIF (nf.EQ.3) THEN
2082  WRITE(nlog2,643,err=113)(vectnb(6*ni+k),k=0,2)
2083  ELSEIF (nf.EQ.4) THEN
2084  WRITE(nlog2,644,err=113)(vectnb(6*ni+k),k=0,3)
2085  ELSEIF (nf.EQ.5) THEN
2086  WRITE(nlog2,645,err=113)(vectnb(6*ni+k),k=0,4)
2087  ENDIF
2088  ENDIF
2089  ENDDO ! FIN BOUCLE EN J
2090  WRITE(nlog2,57,err=113)
2091  CLOSE(nlog2)
2092  ENDDO ! BOUCLE SUR LES SOUS-DOMAINES
2093  CALL system_clock(count=temps_sc(9),count_rate=parsec)
2094  WRITE(lu,*)' REMPLISSAGE DES FICHIERS UNV ET LOG',
2095  & (1.0*(temps_sc(9)-temps_sc(8)))/(1.0*parsec),' SECONDS'
2096 !----------------------------------------------------------------------
2097 ! 6. AFFICHAGES DANS PARTEL.LOG ET TEST DE COMPLETUDE DU PARTITIONNEMENT
2098 !------------
2099 
2100  WRITE(lu,*)' '
2101  compt1=0
2102  compt2=0
2103  compt3=0
2104  DO idd=1,nparts
2105  WRITE(lu,86)idd,nelemsd(idd),npointsd(idd),npointisd(idd)
2106  compt3=compt3+npointisd(idd)
2107  compt2=compt2+npointsd(idd)
2108  compt1=compt1+nelemsd(idd)
2109  ENDDO
2110  WRITE(lu,*)' ------------------------------------'
2111  WRITE(lu,87)compt1,compt2,compt3
2112  WRITE(lu,88)compt1/nparts,compt2/nparts,compt3/nparts
2113  WRITE(lu,*)' '
2114  WRITE(lu,83)(1.0*(temps_sc(9)-temps_sc(1)))/(1.0*parsec)
2115  WRITE(lu,*)' ENDING METIS MESH PARTITIONING--------------------+'
2116  WRITE(lu,*)' '
2117  WRITE(lu,*)' WRITING GEOMETRY FILE FOR EACH PROCESSOR'
2118  WRITE(lu,*)' WRITING LOG FILE FOR EACH PROCESSOR'
2119 
2120 !----------------------------------------------------------------------
2121 ! 7. DIVERS
2122 !---------------
2123 
2124 ! 7.A FORMAT DU LOG
2125 !---------------
2126  50 FORMAT(a80) ! LES AUTRES LIGNES
2127 ! 1234567890123456789012345678901234567890123456789
2128  51 FORMAT(' TOTAL NO. OF NODES : ',i10)
2129  52 FORMAT(' TOTAL NO. OF ELEMENTS : ',i10)
2130  523 FORMAT(' TOTAL NO. OF USER-FLUX : ',i10)
2131  53 FORMAT(' TOTAL NO. OF FAMILIES : ',i10)
2132  531 FORMAT(' TOTAL NUMBER OF EXTERNAL FACES : ',i10)
2133  54 FORMAT(' DEBUT DE KNOLG: ',i10,' ',i10)
2134 
2135  540 FORMAT(6i10) ! LIGNE DE BLOC KNOLG ET PRIORITY
2136  541 FORMAT(i10) ! DERNIERE LIGNE DE BLOC KNOLG
2137  542 FORMAT(2i10) ! DERNIERE LIGNE DE BLOC KNOLG
2138  543 FORMAT(3i10) ! DERNIERE LIGNE DE BLOC KNOLG
2139  544 FORMAT(4i10) ! DERNIERE LIGNE DE BLOC KNOLG
2140  545 FORMAT(5i10) ! DERNIERE LIGNE DE BLOC KNOLG
2141 
2142  641 FORMAT(i9) ! DERNIERE LIGNE DE BLOC NACHB
2143  642 FORMAT(2i9) ! DERNIERE LIGNE DE BLOC NACHB
2144  643 FORMAT(3i9) ! DERNIERE LIGNE DE BLOC NACHB
2145  644 FORMAT(4i9) ! DERNIERE LIGNE DE BLOC NACHB
2146  645 FORMAT(5i9) ! DERNIERE LIGNE DE BLOC NACHB
2147  640 FORMAT(6i9) ! DERNIERE LIGNE DE BLOC NACHB
2148 
2149 
2150 
2151  55 FORMAT(' FIN DE KNOLG: ',i10)
2152  56 FORMAT(' DEBUT DE NACHB: ',i10)
2153  57 FORMAT(' FIN DE NACHB: ')
2154 
2155 ! 7B. FORMAT DU UNV
2156 !---------------
2157  60 FORMAT(a4,a2) ! ' -1'
2158  61 FORMAT(i9) ! LECTURE NSEC
2159  62 FORMAT(a80) ! LECTURE TITRE
2160  63 FORMAT(4i10) ! LIGNE 1 BLOC COORD
2161  64 FORMAT(3d25.16) ! LIGNE 2 BLOC COORD
2162  65 FORMAT(6i10) ! LIGNE 1 BLOC CONNECTIVITE
2163  66 FORMAT(4i10) ! LIGNE 2 BLOC CONNECTIVITE SI TETRA
2164  67 FORMAT(3i10) ! LIGNE 2 BLOC CONNECTIVITE SI TRIANGLE
2165 ! 68 FORMAT(8I10) ! BLOC FANTOCHE POUR MARQUER LA FIN DU BLOC
2166  ! CONNECTIVITEE
2167 
2168 ! 7.C AFFICHAGES DANS PARTEL.LOG
2169 !---------------
2170  80 FORMAT(' #NUMBER TOTAL OF ELEMENTS: ',i8,
2171  & ' #NODES : ',i8)
2172  81 FORMAT(' #TETRAHEDRONS : ',i8,
2173  & ' #TRIANGLE MESH BORDER : ',i8)
2174  82 FORMAT(' #NPARTS : ',i8)
2175  83 FORMAT(' RUNTIME : ',f10.2,' S')
2176  86 FORMAT(' DOMAIN: ',i3,' #ELEMENTS: ',i8,' #NODES: ',i8,
2177  & ' #INTERFACENODES: ',i8)
2178  87 FORMAT(' TOTAL VALUES OF ELEMENTS: ',i10,' NODES: ',i10,
2179  & ' INTERFACENODES: ',i10)
2180  88 FORMAT(' MEAN VALUES OF ELEMENTS : ',i8,' NODES: ',i8,
2181  & ' INTERFACENODES: ',i8)
2182  89 FORMAT(' INPUT UNV FILE :',a50)
2183 ! 90 FORMAT(' INPUT LOG FILE :',A50)
2184 ! 91 FORMAT(' NUMBER OF PARTITIONS:',I5)
2185  92 FORMAT(' NUMBER OF NODES:',i10)
2186  93 FORMAT(' NUMBER OF ELEMENTS:',i10)
2187  94 FORMAT(' NUMBER OF COLORS:',i5)
2188 
2189 ! 7.D DEALLOCATE
2190 !---------------
2191  DEALLOCATE(x1,y1,z1)
2192  DEALLOCATE(ecolor)
2193  DEALLOCATE(iklestet,iklestri,tettri,tettri2)
2194  DEALLOCATE(epart,npart)
2195  DEALLOCATE(nelemsd,npointsd,npointisd)
2196  DEALLOCATE(nodes1,nodes2,nodes3,nodes4,triunv)
2197  DEALLOCATE(nodes1t,nodes2t,nodes3t)
2198  DEALLOCATE(knolg,nachb,priority,ncolor2)
2199  DEALLOCATE(elegl)
2200  DEALLOCATE(nodegl)
2201  DEALLOCATE(nodelg)
2202  DEALLOCATE(nelem_p)
2203  DEALLOCATE(npoin_p)
2204  RETURN
2205 
2206 ! 7.E MESSAGES D'ERREURS
2207 !---------------
2208 
2209  111 texterror='! UNEXPECTED FILE FORMAT2: '//nameinp//' !'
2210  GOTO 999
2211  112 texterror='! UNEXPECTED FILE FORMAT3: '//nameinp2//' !'
2212  GOTO 999
2213  113 texterror='! UNEXPECTED FILE FORMAT4: '//namelog2//' !'
2214  GOTO 999
2215  120 texterror='! UNEXPECTED EOF WHILE READING: '//namelog//' !'
2216  GOTO 999
2217  130 texterror='! PROBLEM WHILE OPENING: '//namelog//' !'
2218  GOTO 999
2219  131 texterror='! PROBLEM WHILE OPENING: '//nameinp//' !'
2220  GOTO 999
2221  132 texterror='! PROBLEM WHILE OPENING: '//nameinp2//' !'
2222  GOTO 999
2223  133 texterror='! PROBLEM WHILE OPENING: '//namelog2//' !'
2224  GOTO 999
2225  140 texterror='! FILE DOES NOT EXIST: '//nameinp//' !'
2226  GOTO 999
2227  141 texterror='! FILE DOES NOT EXIST: '//namelog//' !'
2228  GOTO 999
2229  144 WRITE(unit=str8,fmt='(I8)')maxlensoft
2230  texterror='! NAME OF INPUT FILE '//nameinp//' IS LONGER THAN '//
2231  & str8(1:3)//' CHARACTERS !'
2232  GOTO 999
2233  145 WRITE(unit=str8,fmt='(I8)')maxlensoft
2234  texterror='! NAME OF INPUT FILE '//namelog//' IS LONGER THAN '//
2235  & str8(1:3)//' CHARACTERS !'
2236  GOTO 999
2237  146 texterror='! PROBLEM WITH CONSTRUCTION OF INVERSE CONNECTIVITY !'
2238  GOTO 999
2239  147 texterror='! PROBLEM WHILE WRITING: '//nameinp2//' !'
2240  GOTO 999
2241  149 texterror='! NO INPUT UNV FILE !'
2242  GOTO 999
2243  151 WRITE(unit=str8,fmt='(I8)')j
2244  WRITE(unit=str26,fmt='(I3,1X,I3,1X,I3,1X,I3,1X,I3,1X,I3)')
2245  & (nachb(k,j),k=2,nbsdomvois-1),idd
2246  texterror='! NODE '//str8//' BELONGS TO DOMAINS '//str26(1:23)
2247  & //' !'
2248  GOTO 999
2249  152 texterror='! PROBLEM WITH CONSTRUCTION OF VECTNB FOR NACHB !'
2250  GOTO 999
2251  154 texterror='! PROBLEM WITH THE PRIORITY OF COLOR NODES !'
2252  GOTO 999
2253 ! END OF FILE AND FORMAT ERRORS :
2254  1100 texterror='ERREUR DE LECTURE DU FICHIER UNV '//
2255  & 'VIA MESH_CONNECTIVITY'
2256  GOTO 999
2257  1200 texterror='ERREUR DE FIN DE LECTURE DU FICHIER UNV '//
2258  & 'VIA MESH_CONNECTIVITY'
2259  GOTO 999
2260 
2261  999 WRITE(lu,*) texterror
2262 !
2263  END SUBROUTINE pares3d
subroutine voisin31(IFABOR, NELEM, NELMAX, IELM, IKLE, SIZIKL, NPOIN, NBOR, NPTFR, LIHBOR, KLOG, INDPU, IKLESTR, NELEB2)
Definition: voisin31.f:8
subroutine elebd31(NELBOR, NULONE, IKLBOR, IFABOR, NBOR, IKLE, NELEM, NELEB, NELMAX, NPOIN, NPTFR, IELM)
Definition: elebd31.f:8
subroutine pares3d
Definition: pares3d.F:4