Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
SourceIMP/utilities/RocfracPrep/read_patran.f90
Go to the documentation of this file.
1 !*********************************************************************
2 !* Illinois Open Source License *
3 !* *
4 !* University of Illinois/NCSA *
5 !* Open Source License *
6 !* *
7 !* Copyright@2008, University of Illinois. All rights reserved. *
8 !* *
9 !* Developed by: *
10 !* *
11 !* Center for Simulation of Advanced Rockets *
12 !* *
13 !* University of Illinois *
14 !* *
15 !* www.csar.uiuc.edu *
16 !* *
17 !* Permission is hereby granted, free of charge, to any person *
18 !* obtaining a copy of this software and associated documentation *
19 !* files (the "Software"), to deal with the Software without *
20 !* restriction, including without limitation the rights to use, *
21 !* copy, modify, merge, publish, distribute, sublicense, and/or *
22 !* sell copies of the Software, and to permit persons to whom the *
23 !* Software is furnished to do so, subject to the following *
24 !* conditions: *
25 !* *
26 !* *
27 !* @ Redistributions of source code must retain the above copyright *
28 !* notice, this list of conditions and the following disclaimers. *
29 !* *
30 !* @ Redistributions in binary form must reproduce the above *
31 !* copyright notice, this list of conditions and the following *
32 !* disclaimers in the documentation and/or other materials *
33 !* provided with the distribution. *
34 !* *
35 !* @ Neither the names of the Center for Simulation of Advanced *
36 !* Rockets, the University of Illinois, nor the names of its *
37 !* contributors may be used to endorse or promote products derived *
38 !* from this Software without specific prior written permission. *
39 !* *
40 !* THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, *
41 !* EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES *
42 !* OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND *
43 !* NONINFRINGEMENT. IN NO EVENT SHALL THE CONTRIBUTORS OR *
44 !* COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER *
45 !* LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, *
46 !* ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE *
47 !* USE OR OTHER DEALINGS WITH THE SOFTWARE. *
48 !*********************************************************************
49 !* Please acknowledge The University of Illinois Center for *
50 !* Simulation of Advanced Rockets in works and publications *
51 !* resulting from this software or its derivatives. *
52 !*********************************************************************
53 SUBROUTINE read_patran(numvertx2d,dhmin,nprocs)
54 
55  USE linked_list
56  USE meshdata
57  USE linked_list2
58 
59  IMPLICIT NONE
60 
61  INTEGER :: numvertx2d
62  REAL*8 :: dhmin
63 
64 ! Local
65  REAL*8 :: xx,yy,zz,size1,size2,size3,size4,size5,size6
66  REAL*8 :: dt_courant
67  REAL*8, DIMENSION(1:6) :: disflagvalue, tmpval
68  INTEGER :: numdisflag
69 
70  INTEGER :: id
71  INTEGER :: i, j
72  INTEGER :: itype
73  INTEGER :: n1,n2,n3,n4,n5,n6,n7,n8
74  INTEGER :: iface
75  REAL*8 :: value
76  REAL*8 :: press
77 !-- Tempory holding array of partitioned metis arrays
78  INTEGER, ALLOCATABLE, DIMENSION(:) :: npart
79  INTEGER :: edgecut
80  INTEGER :: nprocs
81  INTEGER :: ip
82 !
83  INTEGER :: iaux,iaux1,iaux2
84  INTEGER, ALLOCATABLE, DIMENSION(:) :: imin,nmin,imax,nmax,ninc
85  INTEGER, ALLOCATABLE, DIMENSION(:) :: imap
86  INTEGER, ALLOCATABLE, DIMENSION(:,:) :: lm4node
87  INTEGER :: nnd1, nnd2, nnd3, nnd4, nnd5, nnd6, nnd7, nnd8, nnd9, nnd10
88  INTEGER :: numnp4node
89  integer, dimension(1:6) :: disflag
90  INTEGER :: ifaceid, itmp
91 
92 
93  INTEGER :: procid
94 
95 
96  LOGICAL :: hex8, tet4, tet10, fileexist
97  LOGICAL :: fileexist1, fileexist2, fileexist3, fileexist4, fileexist5, fileexist6
98  integer :: numgeombodies,iunit,ibody
99 
100  INTEGER, POINTER, DIMENSION(:) :: numnp_glb,numel_glb,nummat_glb
101 
102  INTEGER :: ielemtest, icnt, inodetest, minel, maxel,maxnd, minnd, iel10
103 
104  INTEGER, POINTER, DIMENSION(:) :: epart_loc
105 
106  INTEGER, POINTER, dimension(:) :: procdist
107 
108  TYPE(link_ptr_type) :: link
109 
110  INTEGER :: accumel, accumnd, accumnd_old, accumel_old
111 
112  INTEGER :: ierror, iv, kc
113 
114  size1 = 0.d0
115  size2 = 0.d0
116  size3 = 0.d0
117  size4 = 0.d0
118  size5 = 0.d0
119  size6 = 0.d0
120 
121  ierror = 0
122 
123  hex8 = .false.
124  tet4 = .false.
125  tet10 = .false.
126  overlaymesh = .false.
127  interactmesh = .false.
128 
129  print*,'MESH OPTION:'
130  print*,' READING PATRAN MESH'
131  print*,' '
132 
133  itmp = 0
134  maxnumbc_str = 0
135  maxnumbc_mm = 0
136  maxnumbc_th = 0
137 
138  numndhistory = 0
139 
140  ALLOCATE(numelhex2d(1:5,1:nprocs))
141  ALLOCATE(numeltet2d(1:5,1:nprocs))
142 
143  numelhex2d(:,:) = 0
144  numeltet2d(:,:) = 0
145 
146  fileexist1 = .false.
147  fileexist2 = .false.
148  fileexist3 = .false.
149  fileexist4 = .false.
150  fileexist5 = .false.
151  fileexist6 = .false.
152 
153  INQUIRE (file='Modin/'//prefx(1:prefx_lngth)//'.1.out', exist=fileexist1)
154  INQUIRE (file='Modin/'//prefx(1:prefx_lngth)//'.2.out', exist=fileexist2)
155  INQUIRE (file='Modin/'//prefx(1:prefx_lngth)//'.3.out', exist=fileexist3)
156  INQUIRE (file='Modin/'//prefx(1:prefx_lngth)//'.4.out', exist=fileexist4)
157  INQUIRE (file='Modin/'//prefx(1:prefx_lngth)//'.5.out', exist=fileexist5)
158  INQUIRE (file='Modin/'//prefx(1:prefx_lngth)//'.6.out', exist=fileexist6)
159 
160  numgeombodies = 1
161 
162  IF(fileexist6)THEN
163  numgeombodies = 6
164  OPEN(106,file='Modin/'//prefx(1:prefx_lngth)//'.6.out')
165  OPEN(105,file='Modin/'//prefx(1:prefx_lngth)//'.5.out')
166  OPEN(104,file='Modin/'//prefx(1:prefx_lngth)//'.4.out')
167  OPEN(103,file='Modin/'//prefx(1:prefx_lngth)//'.3.out')
168  OPEN(102,file='Modin/'//prefx(1:prefx_lngth)//'.2.out')
169  OPEN(101,file='Modin/'//prefx(1:prefx_lngth)//'.1.out')
170  ELSE IF (fileexist5)THEN
171  numgeombodies = 5
172  OPEN(105,file='Modin/'//prefx(1:prefx_lngth)//'.5.out')
173  OPEN(104,file='Modin/'//prefx(1:prefx_lngth)//'.4.out')
174  OPEN(103,file='Modin/'//prefx(1:prefx_lngth)//'.3.out')
175  OPEN(102,file='Modin/'//prefx(1:prefx_lngth)//'.2.out')
176  OPEN(101,file='Modin/'//prefx(1:prefx_lngth)//'.1.out')
177  ELSE IF (fileexist4)THEN
178  numgeombodies = 4
179  OPEN(104,file='Modin/'//prefx(1:prefx_lngth)//'.4.out')
180  OPEN(103,file='Modin/'//prefx(1:prefx_lngth)//'.3.out')
181  OPEN(102,file='Modin/'//prefx(1:prefx_lngth)//'.2.out')
182  OPEN(101,file='Modin/'//prefx(1:prefx_lngth)//'.1.out')
183  ELSE IF (fileexist3)THEN
184  numgeombodies = 3
185  OPEN(103,file='Modin/'//prefx(1:prefx_lngth)//'.3.out')
186  OPEN(102,file='Modin/'//prefx(1:prefx_lngth)//'.2.out')
187  OPEN(101,file='Modin/'//prefx(1:prefx_lngth)//'.1.out')
188  ELSE IF (fileexist2) THEN
189  numgeombodies = 2
190  OPEN(102,file='Modin/'//prefx(1:prefx_lngth)//'.2.out')
191  OPEN(101,file='Modin/'//prefx(1:prefx_lngth)//'.1.out')
192  ELSE IF (fileexist1) THEN
193  OPEN(101,file='Modin/'//prefx(1:prefx_lngth)//'.1.out')
194  ELSE
195  fileexist = .false.
196 
197  INQUIRE (file='Modin/'//prefx(1:prefx_lngth)//'.pat', exist=fileexist)
198 
199  IF (.NOT. fileexist)THEN
200  print*,'File '//'Modin/'//prefx(1:prefx_lngth)//'.pat NOT found'
201  print*,'...Trying '//'Modin/'//prefx(1:prefx_lngth)//'.out'
202 
203  fileexist = .false.
204  INQUIRE (file='Modin/'//prefx(1:prefx_lngth)//'.out', exist=fileexist)
205 
206  IF (.NOT. fileexist)THEN
207  print*,'ERROR: File '//'Modin/'//prefx(1:prefx_lngth)//'.out NOT found'
208  print*,'ERROR: No Patran Neutral file Found, stopping'
209  stop
210  ELSE
211  print*,'File '//'Modin/'//prefx(1:prefx_lngth)//'.out found'
212 
213  OPEN(101,file='Modin/'//prefx(1:prefx_lngth)//'.out', form='formatted')
214  ENDIF
215  ELSE
216  OPEN(101,file='Modin/'//prefx(1:prefx_lngth)//'.pat', form='formatted')
217  ENDIF
218 
219  ENDIF
220 
221  print*,'No. of Mesh Segment Files Found =', numgeombodies
222 
223  ALLOCATE(numnp_glb(1:numgeombodies),numel_glb(1:numgeombodies),nummat_glb(1:numgeombodies))
224 
225  DO ibody = 1, numgeombodies
226  iunit = 100 + ibody
227 
228 
229 !
230 ! - Packet Type 25: Title Card, Format(i2,8i8)
231  READ(iunit,'()')
232  READ(iunit,'()')
233 
234 ! - Packet Type 26: Summary Data, Format(i2,8i8)
235 ! - 26 ID IV KC N1 N2 N3 N4 N5
236 ! - N1 = number of nodes
237 ! - N2 = number of elements
238 ! - N3 = number of materials
239 ! - N4 = number of Element Properties
240 ! - N5 = number of Coordinate Frames
241 
242  READ(iunit,*) i,i,i,i,numnp_glb(ibody),numel_glb(ibody),nummat_glb(ibody)
243 
244  IF(nummat_glb(ibody).LE.0) nummat_glb(ibody) = 1
245 
246 ! - part of title card
247  READ(iunit,'()')
248 
249  ENDDO
250 
251  numnp_prmry = sum(numnp_glb)
252  numelv_prmry = sum(numel_glb)
253  nummat = maxval(nummat_glb)
254 
255  print*,' No. of elements = ', numelv_prmry
256  print*,' No. of nodes = ', numnp_prmry
257  print*,' No. of Materials = ', nummat
258  print*,' No. of vertex=',numvertx
259 
260  inodetest = 0
261  ielemtest = 0
262 
263  iaux = 0
264 
265  ALLOCATE(epart(1:numelv_prmry))
266 
267  accumnd = 0
268 
269  accumel = 0
270 
271  DO ibody = 1, numgeombodies
272  iunit = 100 + ibody
273  numnp4node = 0
274 !
275 ! -- Read Nodal coordinates
276 
277  IF(ibody.EQ.1)THEN
278  ALLOCATE(coor(1:3,1:numnp_prmry),lmelv_prmry(1:numvertx,1:numelv_prmry))
279 
280 !!$ IF(numvertx.NE.4)THEN
281 !!$ ALLOCATE(lm4node( 1:4,1:numelv_prmry))
282 !!$ ALLOCATE(imap(1:numnp_prmry))
283 !!$ imap(:) = 0
284 !!$ numnp4node = 0
285 !!$ ELSE
286 !!$ numnp4node = numnp_prmry
287 !!$ ENDIF
288  ENDIF
289 
290  iel10 = 0
291 
292 ! ConvertUnit = .0254d0 ! for titan, units in engles
293 
294  DO i = 1, numnp_glb(ibody)
295 
296 ! - Type 01: Node Data, format(i2,8i8)
297 ! 1 ID IV KC
298 ! ID = node id
299 ! IV = 0 n/a
300 ! KC = 2
301  READ(iunit,*) j,id
302  id = id + accumnd
303 ! - format(3e16.9)
304 ! - caresian coordinate of Nodes x, y, z
305  READ(iunit,'(3e16.9)') coor(1:3,id)
306  coor(1:3,id) = coor(1:3,id)*convertunit
307 ! - ICF GTYPE NDF CONFIG CID PSPC
308  READ(iunit,'()')
309 
310  inodetest = inodetest + 1
311  ENDDO
312  accumnd_old = accumnd
313  accumnd = accumnd + numnp_glb(ibody)
314 !
315 ! -- Read element connectivity array
316  IF(ibody.EQ.1)THEN
317  ALLOCATE(matid(1:numelv_prmry))
318  ALLOCATE(eltypeid(1:numelv_prmry))
319  ENDIF
320 ! flages for element type id
321 !
322 ! 4 = 4 node tetrahedral
323 ! 8 = 8 node hexahedral
324 ! 10 = 10 node tetrahedral
325 ! 12 = 12 node wedge*
326 ! 15 = 15 node wedge*
327 ! 20 = 20 node hexahdral*
328 !
329 ! * Not yet implemented
330 
331  minel = 1e9
332  maxel = 0
333 
334  element: DO i = 1, numel_glb(ibody)
335 
336 ! - Type 02: Element Data, format(i2,8i8)
337 ! 2 ID IV KC N1
338 ! ID = element ID
339 ! IV = Shape (5 = tet)
340  READ(iunit,*) j, id
341  ielemtest = ielemtest + 1
342  accumel_old = accumel
343  id = accumel + id
344 
345 
346 ! NODES CONFIG CID CEID
347  READ(iunit,'(3i8)') numvertx,j, matid(id)
348  IF(numvertx.NE.4)THEN
349  IF(numvertx.NE.10)THEN
350  IF(numvertx.NE.8)THEN
351  print*,'ERROR: FOUND UNSUPPORTED ELEMENT TYPE'
352  print*,'Found element with ', numvertx,' nodes '
353  print*,'Stopping'
354  stop
355  ENDIF
356  ENDIF
357  ENDIF
358 
359  eltypeid(id) = numvertx
360 
361  IF(matid(id).LE.0) matid(id) = 1
362 ! LNODES
363 ! LNODES = Element corner nodes followed by additional nodes
364  READ(iunit,*) lmelv_prmry(1:numvertx,id)
365  lmelv_prmry(1:numvertx,id) = lmelv_prmry(1:numvertx,id) + accumnd_old
366 
367  iel10 = iel10 + 1
368  IF(numvertx.EQ.10.AND.i.EQ.1.AND.nprocs.GT.1)THEN
369  ALLOCATE(lm4node( 1:4,1:numel_glb(ibody)))
370  IF(ibody.EQ.1) ALLOCATE(imap(1:numnp_prmry))
371  imap(:) = 0
372  print*,'Re-setting imap'
373  ENDIF
374  IF(numvertx.EQ.4.AND.i.EQ.1.AND.nprocs.GT.1)THEN
375  ALLOCATE(lm4node( 1:4,1:numel_glb(ibody)))
376  IF(ibody.EQ.1) ALLOCATE(imap(1:numnp_prmry))
377  imap(:) = 0
378  print*,'Re-setting imap'
379  ENDIF
380  IF(numvertx.EQ.10.AND.nprocs.GT.1)THEN
381  DO j = 1, 4
382  IF(imap(lmelv_prmry(j,id)).EQ.0)THEN
383  numnp4node = numnp4node + 1
384  imap(lmelv_prmry(j,id)) = numnp4node
385  ENDIF
386  iaux = iaux + 1
387  lm4node(j,iel10) = imap(lmelv_prmry(j,id))
388  ENDDO
389 
390  ELSE IF(numvertx.EQ.4.AND.nprocs.GT.1)THEN
391  DO j = 1, 4
392  IF(imap(lmelv_prmry(j,id)).EQ.0)THEN
393  numnp4node = numnp4node + 1
394  imap(lmelv_prmry(j,id)) = numnp4node
395  ENDIF
396  iaux = iaux + 1
397  lm4node(j,iel10) = imap(lmelv_prmry(j,id))
398  ENDDO
399  ENDIF
400  maxel = sum(numel_glb(1:ibody)) ! SUM( MAX(MaxEl,id)
401  minel = accumel_old + 1 !MIN(,id)
402 
403 
404 !
405 ! -- Find the size of the smallest element
406 !
407  xx = coor(1,lmelv_prmry(1,id)) - coor(1,lmelv_prmry(2,id))
408  yy = coor(2,lmelv_prmry(1,id)) - coor(2,lmelv_prmry(2,id))
409  zz = coor(3,lmelv_prmry(1,id)) - coor(3,lmelv_prmry(2,id))
410  size1 = sqrt(xx*xx+yy*yy+zz*zz)
411  xx = coor(1,lmelv_prmry(2,id)) - coor(1,lmelv_prmry(3,id))
412  yy = coor(2,lmelv_prmry(2,id)) - coor(2,lmelv_prmry(3,id))
413  zz = coor(3,lmelv_prmry(2,id)) - coor(3,lmelv_prmry(3,id))
414  size2 = sqrt(xx*xx+yy*yy+zz*zz)
415  xx = coor(1,lmelv_prmry(3,id)) - coor(1,lmelv_prmry(1,id))
416  yy = coor(2,lmelv_prmry(3,id)) - coor(2,lmelv_prmry(1,id))
417  zz = coor(3,lmelv_prmry(3,id)) - coor(3,lmelv_prmry(1,id))
418  size3 = sqrt(xx*xx+yy*yy+zz*zz)
419  xx = coor(1,lmelv_prmry(4,id)) - coor(1,lmelv_prmry(1,id))
420  yy = coor(2,lmelv_prmry(4,id)) - coor(2,lmelv_prmry(1,id))
421  zz = coor(3,lmelv_prmry(4,id)) - coor(3,lmelv_prmry(1,id))
422  size4 = sqrt(xx*xx+yy*yy+zz*zz)
423  xx = coor(1,lmelv_prmry(4,id)) - coor(1,lmelv_prmry(2,id))
424  yy = coor(2,lmelv_prmry(4,id)) - coor(2,lmelv_prmry(2,id))
425  zz = coor(3,lmelv_prmry(4,id)) - coor(3,lmelv_prmry(2,id))
426  size5 = sqrt(xx*xx+yy*yy+zz*zz)
427  xx = coor(1,lmelv_prmry(4,id)) - coor(1,lmelv_prmry(3,id))
428  yy = coor(2,lmelv_prmry(4,id)) - coor(2,lmelv_prmry(3,id))
429  zz = coor(3,lmelv_prmry(4,id)) - coor(3,lmelv_prmry(3,id))
430  size6 = sqrt(xx*xx+yy*yy+zz*zz)
431  dhmin = min(size1,size2,size3,size4,size5,size6,dhmin)
432  ENDDO element
433 
434  print*,'Done reading Element'
435  accumel = accumel + numel_glb(ibody)
436 
437  IF(nprocs.GT.1)THEN
438  IF(numvertx.EQ.10)THEN
439  ALLOCATE(npart(1:numnp4node),stat=ierror)
440  IF(ierror .NE.0)THEN
441  ! Space for npart could not be allcoated
442  print*,'Program could not allocate space for npart'
443  print*,'Dimensions 1 to ', numnp4node
444  stop
445  ENDIF
446 
447 !!$ IF(ibody.EQ.4) THEN
448 !!$ WRITE(44,*) NumEL_glb(ibody),numnp4node,nprocs
449 !!$ MaxNd = MAXVAL(lm4node)
450 !!$ MinNd = MINVAL(lm4node)
451 !!$ PRINT*,'MaxNd =', MaxNd
452 !!$ PRINT*,'MinND =', MinNd
453 !!$ PRINT*,'MinNDlocation',MINLOC(lm4node)
454 !!$ PRINT*,'Maximum renumbered node =', numnp4node
455 !!$ PRINT*,'Number of elements in section',NumEL_glb(ibody)
456 !!$ PRINT*,'MinEl:MaxEl',MinEl,MaxEl
457 !!$ DO i = 1, NumEL_glb(ibody)
458 !!$ WRITE(44,*) lm4node(:,i)
459 !!$ IF(ANY(lm4node(:,i) .EQ. 0)) THEN
460 !!$ PRINT*,'zero connectivity'
461 !!$ STOP
462 !!$ ENDIF
463 !!$ ENDDO
464 !!$ CLOSE(44)
465 !!$ ENDIF
466 ! PRINT*,'METIS for elements ',MinEL,':',MaxEl
467 
468 
469 !!$ CALL METIS_PartMeshDual(NumEL_glb(ibody),numnp4node, &
470 !!$ lm4node,2,1,nprocs,edgecut,epart(MinEl:MaxEl),npart)
471 
472  print*,'METIS for elements ',minel,':',maxel,numel_glb(ibody)
473 
474  ALLOCATE(epart_loc(1:numel_glb(ibody)),stat=ierror)
475 
476  IF(ierror .NE.0)THEN
477  ! Space for npart could not be allcoated
478  print*,'Program could not allocate space for epart_loc'
479  print*,'Dimensions 1 to ', numel_glb(ibody)
480  stop
481  ENDIF
482 
483  print*,'CALLING METIS'
484 
485  CALL metis_partmeshnodal(numel_glb(ibody),numnp4node, &
486  lm4node,2,1,nprocs,edgecut,epart_loc,npart)
487 
488  epart(minel:maxel) = epart_loc(1:numel_glb(ibody))
489  print*,minel,maxel,1,numel_glb(ibody)
490  deallocate(epart_loc)
491 
492  print*,'Finished Metis'
493  DEALLOCATE(npart)
494  DEALLOCATE(lm4node)
495 
496  ELSE IF(numvertx.EQ.4)THEN
497  ALLOCATE(npart(1:numnp_glb(ibody)),stat=ierror)
498  IF(ierror .NE.0)THEN
499  ! Space for npart could not be allcoated
500  print*,'Program could not allocate space for npart'
501  print*,'Dimensions 1 to ', numnp4node
502  stop
503  ENDIF
504 
505 !!$ CALL METIS_PartMeshDual(NumEL_glb(ibody),numnp4node, &
506 !!$ lm4node,2,1,nprocs,edgecut,epart(MinEl:MaxEl),npart)
507 
508  print*,'METIS for elements ',minel,':',maxel,numel_glb(ibody)
509 
510  ALLOCATE(epart_loc(1:numel_glb(ibody)),stat=ierror)
511 
512  IF(ierror .NE.0)THEN
513  ! Space for npart could not be allcoated
514  print*,'Program could not allocate space for epart_loc'
515  print*,'Dimensions 1 to ', numel_glb(ibody)
516  stop
517  ENDIF
518 
519  print*,'CALLING METIS',numel_glb(ibody),numnp_glb(ibody)
520 
521  CALL metis_partmeshnodal(numel_glb(ibody),numnp_glb(ibody), &
522  lm4node,2,1,nprocs,edgecut,epart_loc,npart)
523 
524  epart(minel:maxel) = epart_loc(1:numel_glb(ibody))
525  print*,minel,maxel,1,numel_glb(ibody)
526  deallocate(epart_loc)
527 
528  print*,'Finished Metis'
529  DEALLOCATE(npart)
530  DEALLOCATE(lm4node)
531 
532 
533  ENDIF
534  ELSE
535  epart = 1
536  ENDIF
537  ENDDO
538 
539 !!$
540 !!$ OPEN(45,file='LoadBalanceStats.out')
541 !!$ ALLOCATE(ProcDist(1:nprocs))
542 !!$ ProcDist(:) = 0
543 !!$ DO i = 1, numelv_prmry
544 !!$ IF(epart(i).EQ.0)THEN
545 !!$ print*,'epart = 0',i
546 !!$ stop
547 !!$ endif
548 !!$ ProcDist(epart(i)) = ProcDist(epart(i)) + 1
549 !!$ ENDDO
550 !!$
551 !!$ DO i=1, nprocs
552 !!$ write(45,*) i, ProcDist(i)
553 !!$ enddo
554 !!$ PRINT*,'*** Load Balancing Stats ***'
555 !!$ PRINT*,'Maximum number of elements on a processor =', MAXVAL(ProcDist(:))
556 !!$ PRINT*,'Minimum number of elements on a processor =', MINVAL(ProcDist(:))
557 !!$ close(45)
558 !!$ DEALLOCATE(ProcDist)
559 
560 
561 
562  IF(inodetest.NE.numnp_prmry)THEN
563  print*,'ERROR number of combined nodes not equal total nodes'
564  print*,'Sum of Nodes =', numnp_prmry
565  print*,'Actual Number Read =', inodetest
566  stop
567  ENDIF
568 
569  IF(ielemtest.NE.numelv_prmry)THEN
570  print*,'ERROR number of combined elements not equal total elements'
571  print*,'Sum of Elements =', numelv_prmry
572  print*,'Actual Number Read =', ielemtest
573  stop
574  ENDIF
575 
576 ! -- Partition the finite element mesh
577 
578  print*,' REQUESTED NUMBER OF PARTITIONS =',nprocs
579 
580 !
581 ! Call METIS using the partition the mesh
582 !
583 
584 
585  print*,'CALLING METIS'
586 
587  IF(nprocs.GT.1)THEN
588 !-ALT CALL METIS_PartMeshNodal(numel_old_z,nn,
589 !-ALT$ elmnts,1,1,nprocs,edgecut,epart_p,npart_p)
590  IF(numvertx.EQ.8)THEN
591  ALLOCATE(npart(1:numnp_prmry))
592  CALL metis_partmeshdual(numelv_prmry,numnp_prmry, &
593  lmelv_prmry,3,1,nprocs,edgecut,epart,npart)
594 ! DEALLOCATE(npart)
595  ENDIF
596 !!$ ELSE IF(numvertx.EQ.4)THEN
597 !!$ ALLOCATE(npart(1:numnp_prmry))
598 !!$ CALL METIS_PartMeshDual(numelv_prmry,numnp_prmry, &
599 !!$ lmelv_prmry,2,1,nprocs,edgecut,epart,npart)
600 !!$ DEALLOCATE(npart)
601 !!$ ENDIF
602  ELSE
603  epart(:) = 1
604  ENDIF
605 
606  print*,' Called METIS'
607  OPEN(45,file='LoadBalanceStats.out')
608  ALLOCATE(procdist(1:nprocs))
609  procdist(:) = 0
610  DO i = 1, numelv_prmry
611  IF(epart(i).EQ.0)THEN
612  print*,'epart = 0',i
613  stop
614  endif
615  procdist(epart(i)) = procdist(epart(i)) + 1
616  ENDDO
617 
618  DO i=1, nprocs
619  write(45,*) i, procdist(i)
620  enddo
621  print*,'*** Load Balancing Stats ***'
622  print*,'Maximum number of elements on a processor =', maxval(procdist(:))
623  print*,'Minimum number of elements on a processor =', minval(procdist(:))
624  close(45)
625  DEALLOCATE(procdist)
626 
627 
628 ! create communication file ! Note: should be able to remove restriction
629 
630  ALLOCATE(procndlist(1:numnp_prmry,1:maxnumberofprocstosharenode) )
631  ALLOCATE(numprocpernd(1:numnp_prmry))
632  ALLOCATE(numelperproc(1:nprocs))
633  ALLOCATE(numndperproc(1:nprocs))
634 
635  numndperproc(:) = 0
636  numelperproc(:) = 0
637 
638  print*,'Calling NewCommList'
639  CALL newcommlist(numelv_prmry, numnp_prmry, nprocs, numvertx, &
640  lmelv_prmry, epart,numprocpernd,procndlist,maxnumberofprocstosharenode,numelperproc, &
641  numndperproc)
642 
643  print*,' .. completed'
644 
645 !
646 ! Set up the numbering of nodes for the implicit solver
647 !
648 
649  IF ( imp ) THEN
650 
651  ALLOCATE(startnumnp_loc_implicit(1:nprocs))
652  ALLOCATE(numnp_loc_implicit(1:nprocs))
653  ALLOCATE(mapnodeimp(1:numnp_prmry))
654  ALLOCATE(nodeprocimpglobal(1:numnp_prmry))
655  numnp_loc_implicit(:) = 0
656 
657  DO i = 1, numnp_prmry
658  numnp_loc_implicit(npart(i)) = numnp_loc_implicit(npart(i)) + 1
659  mapnodeimp(i) = numnp_loc_implicit(npart(i))
660  nodeprocimpglobal(i) = npart(i) - 1
661  END DO
662 
663  startnumnp_loc_implicit(1) = 0
664  DO i = 2, nprocs
665  startnumnp_loc_implicit(i) = sum(numnp_loc_implicit(1:i-1))
666  ENDDO
667 
668  DO i = 1, numnp_prmry
669  mapnodeimp(i) = mapnodeimp(i) + startnumnp_loc_implicit(npart(i))
670  END DO
671 
672  END IF
673 
674  IF(nprocs.GT.1)THEN
675  IF(numvertx.EQ.8)THEN
676  DEALLOCATE(npart)
677  END IF
678  END IF
679 
680 
681 ! -- create a view of the partitioned mesh
682 ! -- PMVIS software
683 ! -- Command: pmvis.sgi.bin -n pmvis.nod -c pmvis.ele -p pmvis.part -o 1
684 
685  iopmvis = 0
686 
687  IF(iopmvis.EQ.1)THEN
688 
689  OPEN(13,file='pmvis.nod')
690  DO i = 1, numnp_prmry
691  WRITE(13,'(3e16.9)') coor(1:3,i)
692  END DO
693 
694  CLOSE(13)
695  OPEN(13,file='pmvis.ele')
696  DO i = 1, numelv_prmry
697  WRITE(13,'(4i10)') lmelv_prmry(1:4,i)
698  ENDDO
699  CLOSE(13)
700 
701  OPEN(13,file='pmvis.part')
702  DO i = 1, numelv_prmry
703  WRITE(13,'(1i10)') epart(i)
704  ENDDO
705  CLOSE(13)
706  ENDIF
707 
708 !
709 !
710 ! -- Continue to read in the rest of the input parameters
711 
712  numbc_prmry = 0
713  numbc_prmry_mm = 0
714  numbc_prmry_ht = 0
715  numserbc = 0
716 
717 !!$ ALLOCATE(numel_2d(1:nprocs,1:2))
718 !!$ numel_2d(:,:) = 0
719 
720 
721  ALLOCATE(indsburnflg(1:numnp_prmry))
722 
723  indsburnflg(:) = 0
724 
725 ! Initialize list for boundary conditions
726 
727  nullify(bc_structural_head,bc_structural_tail)
728  nullify(bc_meshmotion_head,bc_meshmotion_tail)
729  nullify(bc_thermal_head,bc_thermal_tail)
730 
731  nullify(surfmesh_tri3_ov1_head,surfmesh_tri3_ov1_tail)
732  nullify(surfmesh_tri3_ov2_head,surfmesh_tri3_ov2_tail)
733  nullify(surfmesh_tri3_s_head,surfmesh_tri3_s_tail)
734  nullify(surfmesh_tri3_sf_head,surfmesh_tri3_sf_tail)
735  nullify(surfmesh_tri6_s_head,surfmesh_tri6_s_tail)
736  nullify(surfmesh_tri6_sf_head,surfmesh_tri6_sf_tail)
737  nullify(surfmesh_hex8_s_head,surfmesh_hex8_s_tail)
738  nullify(surfmesh_hex8_sf_head,surfmesh_hex8_sf_tail)
739 
740 ! number of boundary conditions per processors
741 
742  ALLOCATE(numbc_structural(1:nprocs))
743  ALLOCATE(numbc_meshmotion(1:nprocs))
744  ALLOCATE(numbc_thermal(1:nprocs))
745  ALLOCATE(numbc_flag(1:nprocs))
746  numbc_structural(:) = 0
747  numbc_meshmotion(:) = 0
748  numbc_thermal(:) = 0
749  numbc_flag(:) = 0
750  ALLOCATE(bc_flag(1:3,1:numnp_prmry))
751 
752  bc_flag(:,:) = 0
753 
754 
755  print*,'FINISH READING PATRAN FILE'
756  print*,'Reading Boundary Conditions'
757  accumel = 0
758  accumnd = 0
759  DO ibody = 1, numgeombodies
760  iunit = 100 + ibody
761 
762 
763  DO
764  READ(iunit,*) itype,id,iv,kc
765  IF(itype.EQ.99) EXIT
766 
767  IF(itype.EQ.4)THEN
768  READ(iunit,'()')
769 
770 ! Marks of IntFaceFlag :
771 !
772 ! SolidFluid Interface = 1
773 ! NotASolidFluid Interface = 0
774 
775  ELSE IF(itype.EQ. 3) THEN
776  DO i = 1, 20
777  READ(iunit,'()')
778  ENDDO
779  ELSE IF(itype .EQ. 6)THEN ! pressure loading
780  id = id + accumel
781  ip = epart(id)
782  READ(iunit,'(i1,i1,i1,i6,8i1)')i,i,i,i,n1,n2,n3,n4,n5,n6,n7,n8
783  READ(iunit,*) press
784  intfaceflag = abs(int(press))
785 
786  IF(.NOT.(interactmesh))THEN
787  IF(intfaceflag.GE.0.AND.intfaceflag.LE.2) interactmesh = .true.
788  ENDIF
789 
790 ! tet elements
791 
792  IF(n1.EQ.1.AND.n2.EQ.1.AND.n3.EQ.1 &
793  .AND.n4.EQ.0.AND.n5.EQ.0.AND.n6.EQ.0.AND.n7.EQ.0.AND.n8.EQ.0)THEN
794  IF(eltypeid(id).EQ.4)THEN
795  nnd1 = lmelv_prmry(1,id)
796  nnd2 = lmelv_prmry(3,id)
797  nnd3 = lmelv_prmry(2,id)
798  !itmp = itmp + 1
799  ifaceid = 1
800  ELSE
801  nnd1 = lmelv_prmry(1,id)
802  nnd2 = lmelv_prmry(3,id)
803  nnd3 = lmelv_prmry(2,id)
804  nnd4 = lmelv_prmry(7,id)
805  nnd5 = lmelv_prmry(6,id)
806  nnd6 = lmelv_prmry(5,id)
807  ifaceid = 1
808  ENDIF
809  ELSE IF(n1.EQ.1.AND.n2.EQ.1.AND.n4.EQ.1 &
810  .AND.n3.EQ.0.AND.n5.EQ.0.AND.n6.EQ.0.AND.n7.EQ.0.AND.n8.EQ.0)THEN
811  IF(eltypeid(id).EQ.4)THEN
812  nnd1 = lmelv_prmry(1,id)
813  nnd2 = lmelv_prmry(2,id)
814  nnd3 = lmelv_prmry(4,id)
815  !itmp = itmp + 1
816  ifaceid = 2
817  ELSE
818  nnd1 = lmelv_prmry(1,id)
819  nnd2 = lmelv_prmry(2,id)
820  nnd3 = lmelv_prmry(4,id)
821  nnd4 = lmelv_prmry(5,id)
822  nnd5 = lmelv_prmry(9,id)
823  nnd6 = lmelv_prmry(8,id)
824  ifaceid = 2
825  ENDIF
826  ELSE IF(n2.EQ.1.AND.n3.EQ.1.AND.n4.EQ.1 &
827  .AND.n1.EQ.0.AND.n5.EQ.0.AND.n6.EQ.0.AND.n7.EQ.0.AND.n8.EQ.0)THEN
828  IF(eltypeid(id).EQ.4)THEN
829  nnd1 = lmelv_prmry(2,id)
830  nnd2 = lmelv_prmry(3,id)
831  nnd3 = lmelv_prmry(4,id)
832  !itmp = itmp + 1
833  ifaceid = 3
834  ELSE
835  nnd1 = lmelv_prmry(2,id)
836  nnd2 = lmelv_prmry(3,id)
837  nnd3 = lmelv_prmry(4,id)
838  nnd4 = lmelv_prmry(6,id)
839  nnd5 = lmelv_prmry(10,id)
840  nnd6 = lmelv_prmry(9,id)
841  ifaceid = 3
842  ENDIF
843  ELSE IF(n1.EQ.1.AND.n3.EQ.1.AND.n4.EQ.1 &
844  .AND.n2.EQ.0.AND.n5.EQ.0.AND.n6.EQ.0.AND.n7.EQ.0.AND.n8.EQ.0)THEN
845  IF(eltypeid(id).EQ.4)THEN
846  nnd1 = lmelv_prmry(4,id)
847  nnd2 = lmelv_prmry(3,id)
848  nnd3 = lmelv_prmry(1,id)
849  itmp = itmp + 1
850  ifaceid = 4
851  ELSE
852  nnd1 = lmelv_prmry(4,id)
853  nnd2 = lmelv_prmry(3,id)
854  nnd3 = lmelv_prmry(1,id)
855  nnd4 = lmelv_prmry(10,id)
856  nnd5 = lmelv_prmry(7,id)
857  nnd6 = lmelv_prmry(8,id)
858  ifaceid = 4
859  ENDIF
860 ! Hex elements
861  ELSE IF(n1.EQ.1.AND.n2.EQ.1.AND.n3.EQ.1.AND.n4.EQ.1 &
862  .AND.n5.EQ.0.AND.n6.EQ.0.AND.n7.EQ.0.AND.n8.EQ.0)THEN
863  IF(eltypeid(id).EQ.8)THEN
864  nnd1 = lmelv_prmry(4,id)
865  nnd2 = lmelv_prmry(3,id)
866  nnd3 = lmelv_prmry(2,id)
867  nnd4 = lmelv_prmry(1,id)
868  itmp = itmp + 1
869  ifaceid = 1
870  ENDIF
871  ELSE IF(n1.EQ.0.AND.n2.EQ.0.AND.n3.EQ.0.AND.n4.EQ.0 &
872  .AND.n5.EQ.1.AND.n6.EQ.1.AND.n7.EQ.1.AND.n8.EQ.1)THEN
873  IF(eltypeid(id).EQ.8)THEN
874  nnd1 = lmelv_prmry(5,id)
875  nnd2 = lmelv_prmry(6,id)
876  nnd3 = lmelv_prmry(7,id)
877  nnd4 = lmelv_prmry(8,id)
878  itmp = itmp + 1
879  ifaceid = 2
880  ENDIF
881  ELSE IF(n1.EQ.1.AND.n2.EQ.1.AND.n3.EQ.0.AND.n4.EQ.0 &
882  .AND.n5.EQ.1.AND.n6.EQ.1.AND.n7.EQ.0.AND.n8.EQ.0)THEN
883  IF(eltypeid(id).EQ.8)THEN
884  nnd1 = lmelv_prmry(1,id)
885  nnd2 = lmelv_prmry(2,id)
886  nnd3 = lmelv_prmry(6,id)
887  nnd4 = lmelv_prmry(5,id)
888  itmp = itmp + 1
889  ifaceid = 3
890  ENDIF
891  ELSE IF(n1.EQ.0.AND.n2.EQ.1.AND.n3.EQ.1.AND.n4.EQ.0 &
892  .AND.n5.EQ.0.AND.n6.EQ.1.AND.n7.EQ.1.AND.n8.EQ.0)THEN
893  IF(eltypeid(id).EQ.8)THEN
894  nnd1 = lmelv_prmry(2,id)
895  nnd2 = lmelv_prmry(3,id)
896  nnd3 = lmelv_prmry(7,id)
897  nnd4 = lmelv_prmry(6,id)
898  itmp = itmp + 1
899  ifaceid = 4
900  ENDIF
901  ELSE IF(n1.EQ.0.AND.n2.EQ.0.AND.n3.EQ.1.AND.n4.EQ.1 &
902  .AND.n5.EQ.0.AND.n6.EQ.0.AND.n7.EQ.1.AND.n8.EQ.1)THEN
903  IF(eltypeid(id).EQ.8)THEN
904  nnd1 = lmelv_prmry(3,id)
905  nnd2 = lmelv_prmry(4,id)
906  nnd3 = lmelv_prmry(8,id)
907  nnd4 = lmelv_prmry(7,id)
908  itmp = itmp + 1
909  ifaceid = 5
910  ENDIF
911  ELSE IF(n1.EQ.1.AND.n2.EQ.0.AND.n3.EQ.0.AND.n4.EQ.1 &
912  .AND.n5.EQ.1.AND.n6.EQ.0.AND.n7.EQ.0.AND.n8.EQ.1)THEN
913  IF(eltypeid(id).EQ.8)THEN
914  nnd1 = lmelv_prmry(5,id)
915  nnd2 = lmelv_prmry(8,id)
916  nnd3 = lmelv_prmry(4,id)
917  nnd4 = lmelv_prmry(1,id)
918  itmp = itmp + 1
919  ifaceid = 6
920  ENDIF
921  ELSE
922  print*,'Error in pressure face numbering'
923  stop
924  ENDIF
925 
926 ! BC flags
927 ! ignitable (=1)
928 ! non-ignitable (=0)
929 ! noninteracting (=2)
930 
931 
932 ! 3-node triangle, (via 4-node tetrahedral)
933 
934  IF(eltypeid(id).EQ.4)THEN
935 
936  ALLOCATE(surfmesh_tri3_item)
937 
938  surfmesh_tri3_item%ElemData(1) = nnd1
939  surfmesh_tri3_item%ElemData(2) = nnd2
940  surfmesh_tri3_item%ElemData(3) = nnd3
941  surfmesh_tri3_item%ElemData(4) = id
942 
943  IF(intfaceflag.EQ.1)THEN ! Fluid Structure interaction, with burning
944  CALL add_surfmesh_tri3(surfmesh_tri3_item, surfmesh_tri3_sf_head, surfmesh_tri3_sf_tail)
945  numeltet2d(1,ip) = numeltet2d(1,ip) + 1
946 
947  ELSE IF(intfaceflag.EQ.2)THEN ! fluid Structure interaction, without burning
948  CALL add_surfmesh_tri3(surfmesh_tri3_item, surfmesh_tri3_sf_nonignt_head, surfmesh_tri3_sf_nonignt_tail)
949  numeltet2d(3,ip) = numeltet2d(3,ip) + 1
950  ELSE IF(intfaceflag.EQ.0)THEN ! No fluid structure interaction
951  CALL add_surfmesh_tri3(surfmesh_tri3_item, surfmesh_tri3_s_head, surfmesh_tri3_s_tail)
952  numeltet2d(2,ip) = numeltet2d(2,ip) + 1
953  ELSE IF(intfaceflag.EQ.10)THEN ! No fluid structure interaction for non-matching meshes
954  overlaymesh = .true.
955  surfmesh_tri3_item%ElemData(5) = ifaceid
956  CALL add_surfmesh_tri3(surfmesh_tri3_item, surfmesh_tri3_ov1_head, surfmesh_tri3_ov1_tail)
957  numeltet2d(4,ip) = numeltet2d(4,ip) + 1
958  ELSE IF(intfaceflag.EQ.100)THEN ! No fluid structure interaction for non-matching meshes
959  overlaymesh = .true.
960  surfmesh_tri3_item%ElemData(5) = ifaceid
961  CALL add_surfmesh_tri3(surfmesh_tri3_item, surfmesh_tri3_ov2_head, surfmesh_tri3_ov2_tail)
962  numeltet2d(5,ip) = numeltet2d(5,ip) + 1
963 
964  ENDIF
965  tet4 = .true.
966 
967 ! 6-node triangle, (via 10-node tetrahedral)
968 
969  ELSE IF(eltypeid(id).EQ.10)THEN
970 
971  ALLOCATE(surfmesh_tri6_item)
972 
973  surfmesh_tri6_item%ElemData(1) = nnd1
974  surfmesh_tri6_item%ElemData(2) = nnd2
975  surfmesh_tri6_item%ElemData(3) = nnd3
976  surfmesh_tri6_item%ElemData(4) = nnd4
977  surfmesh_tri6_item%ElemData(5) = nnd5
978  surfmesh_tri6_item%ElemData(6) = nnd6
979  surfmesh_tri6_item%ElemData(7) = id
980 
981  IF(intfaceflag.EQ.1)THEN
982  CALL add_surfmesh_tri6(surfmesh_tri6_item, surfmesh_tri6_sf_head, surfmesh_tri6_sf_tail)
983  numeltet2d(1,ip) = numeltet2d(1,ip) + 1
984  ELSEIF(intfaceflag.EQ.2)THEN
985  CALL add_surfmesh_tri6(surfmesh_tri6_item, surfmesh_tri6_sf_nonignt_head, surfmesh_tri6_sf_nonignt_tail)
986  numeltet2d(3,ip) = numeltet2d(3,ip) + 1
987  ELSE IF(intfaceflag.EQ.0) THEN
988  CALL add_surfmesh_tri6(surfmesh_tri6_item, surfmesh_tri6_s_head, surfmesh_tri6_s_tail)
989  numeltet2d(2,ip) = numeltet2d(2,ip) + 1
990  ENDIF
991  tet10 = .true.
992 
993 ! 4-node quad, (via 8-node hexahedral)
994 
995  ELSE IF(eltypeid(id).EQ.8)THEN
996 
997  ALLOCATE(surfmesh_hex8_item)
998 
999  surfmesh_hex8_item%ElemData(1) = nnd1
1000  surfmesh_hex8_item%ElemData(2) = nnd2
1001  surfmesh_hex8_item%ElemData(3) = nnd3
1002  surfmesh_hex8_item%ElemData(4) = nnd4
1003  surfmesh_hex8_item%ElemData(5) = id
1004  IF(intfaceflag.EQ.1)THEN
1005  CALL add_surfmesh_hex8(surfmesh_hex8_item, surfmesh_hex8_sf_head, surfmesh_hex8_sf_tail)
1006  numelhex2d(1,ip) = numelhex2d(1,ip) + 1
1007  ELSE IF(intfaceflag.EQ.2)THEN
1008  CALL add_surfmesh_hex8(surfmesh_hex8_item, surfmesh_hex8_sf_nonignt_head, surfmesh_hex8_sf_nonignt_tail)
1009  numelhex2d(3,ip) = numelhex2d(3,ip) + 1
1010  ELSE IF(intfaceflag.EQ.0) THEN
1011  CALL add_surfmesh_hex8(surfmesh_hex8_item, surfmesh_hex8_s_head, surfmesh_hex8_s_tail)
1012  numelhex2d(2,ip) = numelhex2d(2,ip) + 1
1013  ELSE
1014  print*,'ERROR, interface flag type not found',intfaceflag
1015  stop
1016  ENDIF
1017 
1018  hex8 = .true.
1019  ENDIF
1020 
1021 
1022 ! numel_2d(ip,IntFaceFlag) = numel_2d(ip,IntFaceFlag) + 1
1023 
1024 
1025  ELSE IF(itype .EQ. 8)THEN ! boundary conditions
1026  id = id + accumnd
1027 !
1028 ! CID = Coordinate frame ID
1029 ! ICOMP = 6 displacement component flags (0 or 1)
1030 ! note the flag is passed to node number one (i.e. fixed displacement)
1031 
1032  READ(iunit,'(I8,6I1)') iaux, disflag(1:6)
1033  numdisflag = sum(disflag(1:6))
1034 
1035 
1036  numserbc = numserbc + numdisflag
1037 
1038  READ(iunit,*) tmpval(1:numdisflag)
1039 
1040  icnt = 0
1041  do i = 1, 6
1042  IF(disflag(i).eq.1)THEN
1043  icnt = icnt + 1
1044  disflagvalue(i) = tmpval(icnt)
1045  endif
1046  enddo
1047 
1048  IF(icnt.NE.numdisflag)THEN
1049  print*,'ERROR in processing boundary conditions'
1050  stop
1051  ENDIF
1052 
1053 ! add to the number of boundary conditions per processor
1054 
1055  IF(disflag(1).eq.1 )THEN ! structural boundary condition
1056  DO i = 1,numprocpernd(id)
1057  procid = procndlist(id,i)
1058  numbc_structural(procid) = numbc_structural(procid) + 1
1059  IF(bc_flag(1,id).EQ.0) numbc_flag(procid) = numbc_flag(procid) + 1
1060  ENDDO
1061 
1062  bc_flag(1,id) = bc_flag(1,id) + int(disflagvalue(1))
1063 
1064  maxnumbc_str = max(maxnumbc_str,int(disflagvalue(1)))
1065  endif
1066 
1067  IF(disflag(2).eq.1 )THEN ! temperature boundary condition
1068  DO i = 1,numprocpernd(id)
1069  procid = procndlist(id,i)
1070  numbc_thermal(procid) = numbc_thermal(procid) + 1
1071  IF(bc_flag(2,id).EQ.0) numbc_flag(procid) = numbc_flag(procid) + 1
1072  enddo
1073 
1074  bc_flag(2,id) = bc_flag(2,id) + 100*int(disflagvalue(2))
1075 
1076  maxnumbc_th = max(maxnumbc_th,int(disflagvalue(2)))
1077  endif
1078 
1079 
1080  IF(disflag(3).eq.1 )THEN ! mesh motion boundary condition
1081  DO i = 1,numprocpernd(id)
1082  procid = procndlist(id,i)
1083  numbc_meshmotion(procid) = numbc_meshmotion(procid) + 1
1084  IF(bc_flag(3,id).EQ.0) numbc_flag(procid) = numbc_flag(procid) + 1
1085  ENDDO
1086  bc_flag(3,id) = bc_flag(3,id) + 10000*int(disflagvalue(3))
1087 
1088  maxnumbc_mm = max(maxnumbc_mm,int(disflagvalue(3)))
1089  endif
1090 
1091 
1092 
1093  ELSE IF(itype .EQ. 7)THEN ! mesh motion boundary conditions
1094  id = id + accumnd
1095 !
1096 ! CID = Coordinate frame ID
1097 ! ICOMP = 6 displacement component flags (0 or 1)
1098 ! note the flag is passed to node number one (i.e. fixed displacement)
1099 
1100  READ(iunit,'()')
1101  READ(iunit,*) value
1102 
1103 ! ALLOCATE(BC_meshmotion_item)
1104 
1105 ! BC_meshmotion_item%BC_nodeGlb = id
1106 ! BC_meshmotion_item%BC_flagGlb = INT(value)
1107 
1108 ! CALL add_BC(BC_meshmotion_item, BC_meshmotion_head, BC_meshmotion_tail)
1109 
1110 
1111 ! add to the number of boundary conditions per processor
1112 
1113  DO i = 1,numprocpernd(id)
1114  procid = procndlist(id,i)
1115  numbc_meshmotion(procid) = numbc_meshmotion(procid) + 1
1116  IF(bc_flag(3,id).EQ.0) numbc_flag(procid) = numbc_flag(procid) + 1
1117  ENDDO
1118 
1119  bc_flag(3,id) = bc_flag(3,id) + 10000*int(value)
1120 
1121  maxnumbc_mm = max(maxnumbc_mm,int(value))
1122 
1123  ELSE IF(itype .EQ. 5)THEN ! coordinate frames
1124  READ(iunit,'()')
1125  READ(iunit,'()')
1126  READ(iunit,'()')
1127  READ(iunit,'()')
1128  ELSE IF(itype .EQ. 10)THEN ! nodal temperature (corners)
1129  id = id + accumnd
1130  READ(iunit,*) value
1131 
1132 ! ALLOCATE(BC_thermal_item)
1133 
1134 ! BC_thermal_item%BC_nodeGlb = id
1135 ! BC_thermal_item%BC_flagGlb = INT(value)
1136 
1137 ! CALL add_BC(BC_thermal_item, BC_thermal_head, BC_thermal_tail)
1138 
1139 ! add to the number of boundary conditions per processor
1140 
1141  DO i = 1,numprocpernd(id)
1142  procid = procndlist(id,i)
1143  numbc_thermal(procid) = numbc_thermal(procid) + 1
1144  IF(bc_flag(2,id).EQ.0) numbc_flag(procid) = numbc_flag(procid) + 1
1145  ENDDO
1146 
1147  bc_flag(2,id) = bc_flag(2,id) + 100*int(value)
1148 
1149  maxnumbc_th = max(maxnumbc_th,int(value))
1150 
1151  ELSE
1152  DO i = 1, kc
1153  READ(iunit,'()')
1154  ENDDO
1155  ENDIF
1156  ENDDO
1157 
1158  accumel = sum(numel_glb(1:ibody))
1159 
1160  accumnd = sum(numnp_glb(1:ibody))
1161  ENDDO
1162 
1163 ! if negative number then higher order element
1164 
1165 
1166  IF(tet4.AND..NOT.(tet10).AND..NOT.(hex8))THEN
1167  meshtype2d = 3
1168  ELSE IF (.NOT.(tet4).AND.tet10.AND..NOT.(hex8))THEN
1169  meshtype2d = 6
1170  ELSE IF (.NOT.(tet4).AND..NOT.(tet10).AND.hex8)THEN
1171  meshtype2d = 4
1172  ELSE
1173  meshtype2d = 5
1174  ENDIF
1175 
1176 ! PRINT*,'Number of Surface elements, Interacting with burning =', NumEltet2D(1,:)
1177 ! PRINT*,'Number of Surface elements, Interacting with Non-burning =', NumEltet2D(3,:)
1178 ! PRINT*,'Number of Surface elements, Non-Interacting =', NumEltet2D(2,:)
1179 
1180  print*,'BOUNDARY CONDITIONS'
1181  print*,' Total of all types = ',numserbc
1182  print*,' Number of structural = ',maxnumbc_str
1183  print*,' Number of meshmotion = ',maxnumbc_mm
1184  print*,' Number of thermal = ',maxnumbc_th
1185 
1186  IF(numvertx.EQ.4)THEN
1187  numvertx2d = 3
1188  ELSE IF(numvertx.EQ.10)THEN
1189  numvertx2d = 6
1190  ELSE IF(numvertx.EQ.8)THEN
1191  numvertx2d = 4
1192  ENDIF
1193 
1194  DO ibody = 1, numgeombodies
1195  iunit = 100 + ibody
1196  CLOSE(iunit)
1197  enddo
1198 
1199  RETURN
1200 
1201 END SUBROUTINE read_patran
1202 
1203 
subroutine add_surfmesh_tri6(new, head, tail)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE kc
size_t value()
Returns an index corresponding to a free vertex.
Tfloat sum() const
Return the sum of all the pixel values in an image.
Definition: CImg.h:13022
subroutine newcommlist(NumEl, NumNP, NumProcs, NumVertex, NodeConn, ElPart, NumProcPerNd, ProcNdList, MaxNumberOfProcsToShareNode, NumElPerProc, NumNdPerProc)
Vector_n max(const Array_n_const &v1, const Array_n_const &v2)
Definition: Vector_n.h:354
double sqrt(double d)
Definition: double.h:73
subroutine add_surfmesh_tri3(new, head, tail)
blockLoc i
Definition: read.cpp:79
subroutine read_patran(numvertx2d, dhmin, nprocs)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE form
virtual std::ostream & print(std::ostream &os) const
Vector_n min(const Array_n_const &v1, const Array_n_const &v2)
Definition: Vector_n.h:346
j indices j
Definition: Indexing.h:6
unsigned long id(const Leda_like_handle &x)
Definition: Handle.h:107
subroutine add_surfmesh_hex8(new, head, tail)