KTH framework for Nek5000 toolboxes; testing version  0.0.1
connect1.f
Go to the documentation of this file.
1 c-----------------------------------------------------------------------
2  subroutine setup_topo
3 C
4 C Parallel compatible routine to find
5 C connectivity of element structure.
6 C
7 C On Processor 0:
8 C
9 C .Verify right-handedness of elements.
10 C .Verify element-to-element reciprocity of BC's
11 C .Verify correlation between E-E BC's and physical coincidence
12 C .Set rotations
13 C .Determine multiplicity
14 C .Set up direct stiffness summation arrays.
15 C
16 C All Processors:
17 C
18 C .Disperse/Receive BC and MULT temporary data read from preprocessor.
19 C
20 C
21  include 'SIZE'
22  include 'TOTAL'
23  include 'NONCON'
24  include 'SCRCT'
25 c
26  COMMON /scruz/ xm3(lx3,ly3,lz3,lelt)
27  $ , ym3(lx3,ly3,lz3,lelt)
28  $ , zm3(lx3,ly3,lz3,lelt)
29 C
30  common /c_is1/ glo_num(1*lx1*ly1*lz1*lelv)
31  integer*8 glo_num
32  common /ivrtx/ vertex((2**ldim)*lelt)
33  integer vertex
34 
35  if(nio.eq.0) write(6,*) 'setup mesh topology'
36 C
37 C Initialize key arrays for Direct Stiffness SUM.
38 C
39  nxl=3
40  nyl=3
41  nzl=1+2*(ldim-2)
42 
43  call initds
44  call dsset (nxl,nyl,nzl)
45  call setedge
46 C
47 C=================================================
48 C Establish (global) domain topology
49 C=================================================
50 C
51 C .Generate topologically correct mesh data.
52 C .Set up element centers, face centers, etc.
53 C .Check right handedness of elements.
54 C .Check element boundary conditions.
55 C .Establish Element-Element rotations
56 C .Construct the element to processor map and
57 
58  call genxyzl
59  call setside
60  call verify
61 
62  CALL setcdof
63  IF (ifaxis ) CALL setrzer
64  IF (ifmvbd ) CALL cbcmesh
65  CALL chkaxcb
66 C
67 C========================================================================
68 C Set up element-processor mapping and establish global numbering
69 C========================================================================
70 C
71  mfield=2
72  if (ifflow) mfield=1
73  if (ifmvbd) mfield=0
74 
75  ncrnr = 2**ldim
76 
77  if (nelgv.eq.nelgt) then
78  call get_vert
79  call setupds(gsh_fld(1),lx1,ly1,lz1,nelv,nelgv,vertex,glo_num)
80  gsh_fld(2)=gsh_fld(1)
81 
82 c call gs_counter(glo_num,gsh_fld(1))
83 
84  else
85 
86 c
87 c For conjugate heat transfer, it is assumed that fluid
88 c elements are listed both globally and locally with lower
89 c element numbers than the solid elements.
90 c
91 c We currently assume that there is at least one fluid elem.
92 c per processor.
93 c
94 
95  call get_vert
96 c call outmati(vertex,4,nelv,'vrtx V')
97  call setupds(gsh_fld(1),lx1,ly1,lz1,nelv,nelgv,vertex,glo_num)
98 
99 c call get_vert (vertex, ncrnr, nelgt, '.mp2') ! LATER !
100 c call outmati(vertex,4,nelt,'vrtx T')
101  call setupds(gsh_fld(2),lx1,ly1,lz1,nelt,nelgt,vertex,glo_num)
102 
103 c
104 c Feb 20, 2012: It appears that we do not need this restriction: (pff)
105 c
106 c check if there is a least one fluid element on each processor
107 c do iel = 1,nelt
108 c ieg = lglel(iel)
109 c if (ieg.le.nelgv) goto 101
110 c enddo
111 c if(nio.eq.0) write(6,*)
112 c & 'ERROR: each domain must contain at least one fluid element!'
113 c call exitt
114 c 101 continue
115 
116  endif
117 
118 
119 c if (ifmvbd) call setup_mesh_dssum ! Set up dssum for mesh
120 
121 C========================================================================
122 C Set up multiplicity and direct stiffness arrays for each IFIELD
123 C========================================================================
124 
125  ntotv = lx1*ly1*lz1*nelv
126  ntott = lx1*ly1*lz1*nelt
127 
128 
129 
130  if (ifflow) then
131  ifield = 1
132  call rone (vmult,ntotv)
133  call dssum (vmult,lx1,ly1,lz1)
134  vmltmax=glmax(vmult,ntotv)
135  ivmltmax=vmltmax
136  if (nio.eq.0) write(6,*) ivmltmax,' max multiplicity'
137  call invcol1 (vmult,ntotv)
138  endif
139  if (ifheat) then
140  ifield = 2
141  call rone (tmult,ntott)
142  call dssum (tmult,lx1,ly1,lz1)
143  call invcol1 (tmult,ntott)
144  endif
145  if (.not.ifflow) call copy(vmult,tmult,ntott)
146  if (ifmvbd) call copy (wmult,vmult,ntott)
147  do ifield=3,nfield ! Additional pass. scalrs.
148  if (nelg(ifield).eq.nelgv) then
149  gsh_fld(ifield) = gsh_fld(1)
150  call copy (tmult(1,1,1,1,ifield-1),vmult,ntotv)
151  else
152  gsh_fld(ifield) = gsh_fld(2)
153  call copy (tmult(1,1,1,1,ifield-1),tmult,ntott)
154  endif
155  enddo
156 
157  ifgsh_fld_same = .true.
158  do ifield=2,nfield
159  if (gsh_fld(ifield).ne.gsh_fld(1)) then
160  ifgsh_fld_same = .false.
161  goto 100
162  endif
163  enddo
164  100 continue
165 
166  if(nio.eq.0) then
167  write(6,*) 'done :: setup mesh topology'
168  write(6,*) ' '
169  endif
170 
171  return
172  end
173 c-----------------------------------------------------------------------
174  subroutine initds
175 C
176 C -- Direct Stiffness Initialization Routine --
177 C
178 C Set up required data for packing data on faces of spectral cubes.
179 C
180  include 'SIZE'
181  include 'TOPOL'
182 C
183 C Nominal ordering for direct stiffness summation of faces
184 C
185  j=0
186  DO 5 idim=1,ldim
187  DO 5 iface=1,2
188  j=j+1
189  nomlis(iface,idim)=j
190  5 CONTINUE
191 C
192 C Assign Ed's numbering scheme to PF's scheme.
193 C
194  eface(1)=4
195  eface(2)=2
196  eface(3)=1
197  eface(4)=3
198  eface(5)=5
199  eface(6)=6
200 C
201 C Assign inverse of Ed's numbering scheme to PF's scheme.
202 C
203  eface1(1)=3
204  eface1(2)=2
205  eface1(3)=4
206  eface1(4)=1
207  eface1(5)=5
208  eface1(6)=6
209 C
210 C Assign group designation to each face to determine ordering of indices.
211 C
212  group(1)=0
213  group(2)=1
214  group(3)=1
215  group(4)=0
216  group(5)=0
217  group(6)=1
218 C
219  RETURN
220  END
221 c-----------------------------------------------------------------------
222  subroutine setedge
223 C
224 C .Initialize EDGE arrays for face and edge specific tasks.
225 C
226 C .NOTE: Sevaral arrays in common are initialized via
227 C BLOCKDATA EDGEC
228 C
229 C Computed arrays:
230 C
231 C IEDGE - Minimal list of wire frame nodes.
232 C Used to search for all physical
233 C coincidences.
234 C
235 C
236 C IEDGEF - .Ordered list of wire frame nodes
237 C associated with faces 1 through 6.
238 C .Each of 4 sides of square frame stored
239 C individually so that rotations are
240 C readily handled.
241 C .Two types of node orderings stored -
242 C (0) is clockwise marching
243 C (1) is counter-clockwise marching
244 C for image face.
245 C
246 C
247 C IFACE - indicates the face number. Two notations
248 C are currently in use:
249 C
250 C i) Preprocessor notation:
251 C
252 C +--------+ ^ S
253 C / /| |
254 C / 3 / | |
255 C 4--> / / | |
256 C +--------+ 2 + +----> R
257 C | | / /
258 C | 6 | / /
259 C | |/ /
260 C +--------+ T
261 C 1
262 C
263 C ii) Symmetric notation:
264 C
265 C +--------+ ^ S
266 C / /| |
267 C / 4 / | |
268 C 1--> / / | |
269 C +--------+ 2 + +----> R
270 C | | / /
271 C | 6 | / /
272 C | |/ /
273 C +--------+ T
274 C 3
275 C
276 C EFACE(IFACE) - Given face number IFACE in symmetric notation,
277 C returns preprocessor notation face number.
278 C
279 C EFACE1(IFACE) - Given face number IFACE in preprocessor notation,
280 C returns symmetric notation face number.
281 C
282 C The following variables all take the symmetric
283 C notation of IFACE as arguments:
284 C
285 C ICFACE(i,IFACE) - Gives the 4 vertices which reside on face IFACE
286 C as depicted below, e.g. ICFACE(i,2)=2,4,6,8.
287 C
288 C 3+-----+4 ^ Y
289 C / 2 /| |
290 C Edge 1 extends / / | |
291 C from vertex 7+-----+8 +2 +----> X
292 C 1 to 2. | 4 | / /
293 C | |/ /
294 C 5+-----+6 Z
295 C 3
296 C
297 C IEDGFC(i,IFACE) - Gives the 4 edges which border the face IFACE
298 C Edge numbering is as follows:
299 C Edge = 1,2,3,4 run in +r direction
300 C Edge = 5,6,7,8 run in +s direction
301 C Edge = 9,10,11,12 run in +t direction
302 C
303 C Ordering of each edge is such that a monotonically
304 C increasing sequence of vertices is associated with
305 C the start point of a corresponding set of
306 C monotonically increasing edge numbers, e.g.,
307 C
308 C ICEDG(i,IEDGE) - Gives 3 variables for determining the stride along
309 C a given edge, IEDGE; i=1 gives the starting vertex
310 C i=2 gives the stopping vertex
311 C i=3 gives the stride size.
312 C
313  include 'SIZE'
314  include 'TOPOL'
315 C
316  COMMON /ctmp0/ itmp(3,3,3)
317  INTEGER ORDER
318 C
319  nxl=3
320  nyl=3
321  nzl=1+2*(ldim-2)
322  nxy =nxl*nyl
323  nxyz =nxl*nyl*nzl
324  nfaces=2*ldim
325 C
326 C----------------------------------------------------------------------
327 C Set up edge arrays (temporary - required only for defining DS)
328 C----------------------------------------------------------------------
329 C
330 C Fill corners - 1 through 8.
331 C
332  i3d=1
333  IF (ldim.EQ.2) i3d=0
334 C
335  i=0
336  DO 10 i3=0,i3d
337  iz=1+(nzl-1)*i3
338  DO 10 i2=0,1
339  iy=1+(nyl-1)*i2
340  DO 10 i1=0,1
341  ix=1+(nxl-1)*i1
342  i=i+1
343  iedge(i)=ix+nxl*(iy-1)+nxy*(iz-1)
344  10 CONTINUE
345 C
346 C Fill X-direction edges.
347 C
348  DO 20 i3=0,i3d
349  iz=1+(nzl-1)*i3
350  DO 20 i2=0,1
351  iy=1+(nyl-1)*i2
352  DO 20 ix=2,nxl-1
353  i=i+1
354  iedge(i)=ix+nxl*(iy-1)+nxy*(iz-1)
355  20 CONTINUE
356 C
357 C Fill Y-direction edges.
358 C
359  DO 30 i3=0,i3d
360  iz=1+(nzl-1)*i3
361  DO 30 i1=0,1
362  ix=1+(nxl-1)*i1
363  DO 30 iy=2,nyl-1
364  i=i+1
365  iedge(i)=ix+nxl*(iy-1)+nxy*(iz-1)
366  30 CONTINUE
367 C
368 C Fill Z-direction edges.
369 C
370  IF (ldim.EQ.3) THEN
371  DO 40 i2=0,1
372  iy=1+(nyl-1)*i2
373  DO 40 i1=0,1
374  ix=1+(nxl-1)*i1
375  DO 40 iz=2,nzl-1
376  i=i+1
377  iedge(i)=ix+nxl*(iy-1)+nxy*(iz-1)
378  40 CONTINUE
379  ENDIF
380 C
381  CALL izero(invedg,27)
382  DO 44 ii=1,20
383  ix=iedge(ii)
384  invedg(ix)=ii
385  44 CONTINUE
386 C
387 C
388 C GENERAL FACE, GENERAL ROTATION EDGE NUMBERS.
389 C
390  IF (ldim.EQ.3) THEN
391 C
392 C Pack 3-D edge numbering:
393 C
394 C Fill temporary array with local index numbers:
395 C
396  DO 50 ix=1,nxyz
397  itmp(ix,1,1)=ix
398  50 CONTINUE
399 C
400 C Two sets are required, the base cube and the image cube
401 C which is being summed with it.
402 C
403  DO 1000 image=0,1
404 C
405 C Pack edges for each face, no rotation.
406 C
407  DO 500 iface=1,nfaces
408  js1 = skpdat(1,iface)
409  jf1 = skpdat(2,iface)
410  jskip1 = skpdat(3,iface)
411  js2 = skpdat(4,iface)
412  jf2 = skpdat(5,iface)
413  jskip2 = skpdat(6,iface)
414 C
415 C Choose proper indexing order according to face type and image.
416 C
417  order = (-1)**(group(iface)+image)
418  IF (order.EQ.1) THEN
419 C
420 C Forward ordering:
421 C
422 C +-------------+ ^ v1
423 C | --------->| | |
424 C | ^ 2 | | +-->
425 C | | | | v2
426 C | |1 3| |
427 C | | 4 V |
428 C | |<--------- |
429 C F-------------I F is fiducial node.
430 C
431 C I is location of fiducial node for
432 C image face.
433 C
434 C Load edge 1:
435 C
436  j=0
437  j2=js2
438  DO 100 j1=js1,jf1-jskip1,jskip1
439  j=j+1
440  iedgef(j,1,iface,image)=itmp(j1,j2,1)
441  100 CONTINUE
442 C
443 C Load edge 2:
444 C
445  j=0
446  j1=jf1
447  DO 200 j2=js2,jf2-jskip2,jskip2
448  j=j+1
449  iedgef(j,2,iface,image)=itmp(j1,j2,1)
450  200 CONTINUE
451 C
452 C
453 C Load edge 3:
454 C
455  j=0
456  j2=jf2
457  DO 300 j1=jf1,js1+jskip1,-jskip1
458  j=j+1
459  iedgef(j,3,iface,image)=itmp(j1,j2,1)
460  300 CONTINUE
461 C
462 C Load edge 4:
463 C
464  j=0
465  j1=js1
466  DO 400 j2=jf2,js2+jskip2,-jskip2
467  j=j+1
468  iedgef(j,4,iface,image)=itmp(j1,j2,1)
469  400 CONTINUE
470 C
471  ELSE
472 C
473 C Reverse ordering:
474 C
475 C +-------------+
476 C | |<--------- | ^ v2
477 C | | 2 ^ | |
478 C | | | | <--+
479 C | |3 1| | v1
480 C | V 4 | |
481 C | --------->| |
482 C I-------------F F is fiducial node.
483 C
484 C I is location of fiducial node for
485 C image face.
486 C
487 C Load edge 1:
488 C
489  j=0
490  j1=js1
491  DO 105 j2=js2,jf2-jskip2,jskip2
492  j=j+1
493  iedgef(j,1,iface,image)=itmp(j1,j2,1)
494  105 CONTINUE
495 C
496 C Load edge 2:
497 C
498  j=0
499  j2=jf2
500  DO 205 j1=js1,jf1-jskip1,jskip1
501  j=j+1
502  iedgef(j,2,iface,image)=itmp(j1,j2,1)
503  205 CONTINUE
504 C
505 C Load edge 3:
506 C
507  j=0
508  j1=jf1
509  DO 305 j2=jf2,js2+jskip2,-jskip2
510  j=j+1
511  iedgef(j,3,iface,image)=itmp(j1,j2,1)
512  305 CONTINUE
513 C
514 C Load edge 4:
515 C
516  j=0
517  j2=js2
518  DO 405 j1=jf1,js1+jskip1,-jskip1
519  j=j+1
520  iedgef(j,4,iface,image)=itmp(j1,j2,1)
521  405 CONTINUE
522  ENDIF
523 C
524  500 CONTINUE
525  1000 CONTINUE
526  ELSE
527 C
528 C Load edge information for 2-D case
529 C
530  iedgef(1,1,1,0) = nxy - nxl + 1
531  iedgef(1,2,1,0) = 1
532  iedgef(1,1,2,0) = nxl
533  iedgef(1,2,2,0) = nxy
534  iedgef(1,1,3,0) = 1
535  iedgef(1,2,3,0) = nxl
536  iedgef(1,1,4,0) = nxy
537  iedgef(1,2,4,0) = nxy - nxl + 1
538 C
539  iedgef(1,1,1,1) = 1
540  iedgef(1,2,1,1) = nxy - nxl + 1
541  iedgef(1,1,2,1) = nxy
542  iedgef(1,2,2,1) = nxl
543  iedgef(1,1,3,1) = nxl
544  iedgef(1,2,3,1) = 1
545  iedgef(1,1,4,1) = nxy - nxl + 1
546  iedgef(1,2,4,1) = nxy
547  ENDIF
548 C
549  RETURN
550  END
551 c-----------------------------------------------------------------------
552  subroutine dsset(nx,ny,nz)
553 C
554 C Set up arrays IXCN,ESKIP,SKPDAT,NEDG,NOFFST for new NX,NY,NZ
555 C
556  include 'SIZE'
557  include 'INPUT'
558  include 'TOPOL'
559  INTEGER NXO,NYO,NZO
560  SAVE nxo,nyo,nzo
561  DATA nxo,nyo,nzo /3*0/
562 C
563 C Check if element surface counters are already set from last call...
564 C
565  IF (nxo.EQ.nx.AND.nyo.EQ.ny.AND.nzo.EQ.nz) RETURN
566 C
567 C else, proceed....
568 C
569  nxo = nx
570  nyo = ny
571  nzo = nz
572 C
573 C Establish corner to elemental node number mappings
574 C
575  ic=0
576  DO 10 icz=0,1
577  DO 10 icy=0,1
578  DO 10 icx=0,1
579 C Supress vectorization to
580 c IF(ICX.EQ.0)DUMMY=0
581 c IF(ICX.EQ.1)DUMMY=1
582 c DUMMY2=DUMMY2+DUMMY
583  ic=ic+1
584  ixcn(ic)= 1 + (nx-1)*icx + nx*(ny-1)*icy + nx*ny*(nz-1)*icz
585  10 CONTINUE
586 C
587 C Assign indices for direct stiffness summation of arbitrary faces.
588 C
589 C
590 C Y-Z Planes (Faces 1 and 2)
591 C
592  skpdat(1,1)=1
593  skpdat(2,1)=nx*(ny-1)+1
594  skpdat(3,1)=nx
595  skpdat(4,1)=1
596  skpdat(5,1)=ny*(nz-1)+1
597  skpdat(6,1)=ny
598 C
599  skpdat(1,2)=1 + (nx-1)
600  skpdat(2,2)=nx*(ny-1)+1 + (nx-1)
601  skpdat(3,2)=nx
602  skpdat(4,2)=1
603  skpdat(5,2)=ny*(nz-1)+1
604  skpdat(6,2)=ny
605 C
606 C X-Z Planes (Faces 3 and 4)
607 C
608  skpdat(1,3)=1
609  skpdat(2,3)=nx
610  skpdat(3,3)=1
611  skpdat(4,3)=1
612  skpdat(5,3)=ny*(nz-1)+1
613  skpdat(6,3)=ny
614 C
615  skpdat(1,4)=1 + nx*(ny-1)
616  skpdat(2,4)=nx + nx*(ny-1)
617  skpdat(3,4)=1
618  skpdat(4,4)=1
619  skpdat(5,4)=ny*(nz-1)+1
620  skpdat(6,4)=ny
621 C
622 C X-Y Planes (Faces 5 and 6)
623 C
624  skpdat(1,5)=1
625  skpdat(2,5)=nx
626  skpdat(3,5)=1
627  skpdat(4,5)=1
628  skpdat(5,5)=ny
629  skpdat(6,5)=1
630 C
631  skpdat(1,6)=1 + nx*ny*(nz-1)
632  skpdat(2,6)=nx + nx*ny*(nz-1)
633  skpdat(3,6)=1
634  skpdat(4,6)=1
635  skpdat(5,6)=ny
636  skpdat(6,6)=1
637 C
638 C Set up skip indices for each of the 12 edges
639 C
640 C Note that NXY = NX*NY even for 2-D since
641 C this branch does not apply to the 2D case anyway.
642 C
643 C ESKIP(*,1) = start location
644 C ESKIP(*,2) = end
645 C ESKIP(*,3) = stride
646 C
647  nxy=nx*ny
648  eskip( 1,1) = ixcn(1) + 1
649  eskip( 1,2) = ixcn(2) - 1
650  eskip( 1,3) = 1
651  eskip( 2,1) = ixcn(3) + 1
652  eskip( 2,2) = ixcn(4) - 1
653  eskip( 2,3) = 1
654  eskip( 3,1) = ixcn(5) + 1
655  eskip( 3,2) = ixcn(6) - 1
656  eskip( 3,3) = 1
657  eskip( 4,1) = ixcn(7) + 1
658  eskip( 4,2) = ixcn(8) - 1
659  eskip( 4,3) = 1
660  eskip( 5,1) = ixcn(1) + nx
661  eskip( 5,2) = ixcn(3) - nx
662  eskip( 5,3) = nx
663  eskip( 6,1) = ixcn(2) + nx
664  eskip( 6,2) = ixcn(4) - nx
665  eskip( 6,3) = nx
666  eskip( 7,1) = ixcn(5) + nx
667  eskip( 7,2) = ixcn(7) - nx
668  eskip( 7,3) = nx
669  eskip( 8,1) = ixcn(6) + nx
670  eskip( 8,2) = ixcn(8) - nx
671  eskip( 8,3) = nx
672  eskip( 9,1) = ixcn(1) + nxy
673  eskip( 9,2) = ixcn(5) - nxy
674  eskip( 9,3) = nxy
675  eskip(10,1) = ixcn(2) + nxy
676  eskip(10,2) = ixcn(6) - nxy
677  eskip(10,3) = nxy
678  eskip(11,1) = ixcn(3) + nxy
679  eskip(11,2) = ixcn(7) - nxy
680  eskip(11,3) = nxy
681  eskip(12,1) = ixcn(4) + nxy
682  eskip(12,2) = ixcn(8) - nxy
683  eskip(12,3) = nxy
684 C
685 C Load reverse direction edge arrays for reverse mappings...
686 C
687  DO 20 ied=1,12
688  iedm=-ied
689  eskip(iedm,1) = eskip(ied,2)
690  eskip(iedm,2) = eskip(ied,1)
691  eskip(iedm,3) = -eskip(ied,3)
692  20 CONTINUE
693 C
694 C Compute offset for global edge vector given current element
695 C dimensions.....
696 C
697 C NGSPED(ITE,ICMP) = number of global (ie, distinct) special edges
698 C of type ITE (1,2, or 3) for field ICMP.
699 C
700 C ITE = 1 implies an "X" edge
701 C ITE = 2 implies an "Y" edge
702 C ITE = 3 implies an "Z" edge
703 C
704 C Set up number of nodes along each of the 3 types of edges
705 C (endpoints excluded).
706 C
707  nedg(1)=nx-2
708  nedg(2)=ny-2
709  nedg(3)=nz-2
710 C
711 C
712  RETURN
713  END
714 c-----------------------------------------------------------------------
715  subroutine genxyzl
716 C
717 C Generate xyz coordinates
718 C
719  include 'SIZE'
720  include 'INPUT'
721  include 'SCRCT'
722  COMMON /ctmp0/ xcb(2,2,2),ycb(2,2,2),zcb(2,2,2),h(3,3,2),indx(8)
723 C
724  nxl=3
725  nyl=3
726  nzl=1+2*(ldim-2)
727  ntot3=nxl*nyl*nzl*nelt
728 C
729 C Preprocessor Corner notation: Symmetric Corner notation:
730 C
731 C 4+-----+3 ^ s 3+-----+4 ^ s
732 C / /| | / /| |
733 C / / | | / / | |
734 C 8+-----+7 +2 +----> r 7+-----+8 +2 +----> r
735 C | | / / | | / /
736 C | |/ / | |/ /
737 C 5+-----+6 t 5+-----+6 t
738 C
739  DO 10 ix=1,nxl
740  h(ix,1,1)=0.5*float(3-ix)
741  h(ix,1,2)=0.5*float(ix-1)
742  10 CONTINUE
743  DO 20 iy=1,nyl
744  h(iy,2,1)=0.5*float(3-iy)
745  h(iy,2,2)=0.5*float(iy-1)
746  20 CONTINUE
747  DO 30 iz=1,nzl
748  h(iz,3,1)=0.5*float(3-iz)
749  h(iz,3,2)=0.5*float(iz-1)
750  30 CONTINUE
751 C
752  indx(1)=1
753  indx(2)=2
754  indx(3)=4
755  indx(4)=3
756  indx(5)=5
757  indx(6)=6
758  indx(7)=8
759  indx(8)=7
760 C
761  CALL rzero(xml,ntot3)
762  CALL rzero(yml,ntot3)
763  CALL rzero(zml,ntot3)
764  CALL rzero(xcb,8)
765  CALL rzero(ycb,8)
766  CALL rzero(zcb,8)
767 C
768  DO 5000 ie=1,nelt
769 C
770  ldim2 = 2**ldim
771  DO 50 ix=1,ldim2
772  i=indx(ix)
773  xcb(ix,1,1)=xc(i,ie)
774  ycb(ix,1,1)=yc(i,ie)
775  zcb(ix,1,1)=zc(i,ie)
776  50 CONTINUE
777 C
778 C Map R-S-T space into physical X-Y-Z space.
779 C
780  DO 100 izt=1,ldim-1
781  DO 100 iyt=1,2
782  DO 100 ixt=1,2
783 C
784  DO 100 iz=1,nzl
785  DO 100 iy=1,nyl
786  DO 100 ix=1,nxl
787  xml(ix,iy,iz,ie)=xml(ix,iy,iz,ie)+
788  $ h(ix,1,ixt)*h(iy,2,iyt)*h(iz,3,izt)*xcb(ixt,iyt,izt)
789  yml(ix,iy,iz,ie)=yml(ix,iy,iz,ie)+
790  $ h(ix,1,ixt)*h(iy,2,iyt)*h(iz,3,izt)*ycb(ixt,iyt,izt)
791  zml(ix,iy,iz,ie)=zml(ix,iy,iz,ie)+
792  $ h(ix,1,ixt)*h(iy,2,iyt)*h(iz,3,izt)*zcb(ixt,iyt,izt)
793  100 CONTINUE
794 C
795  5000 CONTINUE
796  RETURN
797  END
798 c-----------------------------------------------------------------------
799  subroutine verify
800 C
801 C .Verify right-handedness of elements.
802 C .Verify element-to-element reciprocity of BC's
803 C .Verify correlation between E-E BC's and physical coincidence
804 C
805  include 'SIZE'
806  include 'PARALLEL'
807  include 'INPUT'
808  include 'SCRCT'
809 
810  call verrhe
811 
812  return
813  end
814 c-----------------------------------------------------------------------
815  subroutine setside
816  include 'SIZE'
817  include 'INPUT'
818  include 'TOPOL'
819  include 'SCRCT'
820 C
821 C SIDE(i,IFACE,IE) - Physical (xyz) location of element side midpoint.
822 C i=1,2,3 gives x,y,z value, respectively.
823 C i=4 gives average dimension of face for setting
824 C tolerances.
825 C
826  indx(1)=1
827  indx(2)=2
828  indx(3)=4
829  indx(4)=3
830  indx(5)=5
831  indx(6)=6
832  indx(7)=8
833  indx(8)=7
834 C
835 C Flip vertex array structure
836 C
837 c write(6,*) nelv,nelt,if3d
838  call rzero(xyz,24*nelt)
839  if (if3d) then
840  do ie=1,nelt
841  do j=1,8
842  ivtx = indx(j)
843  xyz(1,ivtx,ie) = xc(j,ie)
844  xyz(2,ivtx,ie) = yc(j,ie)
845  xyz(3,ivtx,ie) = zc(j,ie)
846 c write(6,1) ie,j,ivtx,xc(j,ie),yc(j,ie),zc(j,ie),' xcz'
847 c write(6,1) ie,j,ivtx,(xyz(k,ivtx,ie),k=1,3),' vtx'
848 c 1 format(3i5,1p3e12.4,a4)
849  enddo
850  enddo
851  else
852  do ie=1,nelt
853  do j=1,4
854  ivtx = indx(j)
855  xyz(1,ivtx,ie) = xc(j,ie)
856  xyz(2,ivtx,ie) = yc(j,ie)
857  xyz(3,ivtx,ie) = 0.0
858  enddo
859  enddo
860  endif
861 C
862 C Compute location of center and "diameter" of each element side.
863 C
864  nfaces=ldim*2
865  ncrnr =2**(ldim-1)
866  CALL rzero(side,24*nelt)
867  DO 500 icrn=1,ncrnr
868  DO 500 ifac=1,nfaces
869  ivtx = icface(icrn,ifac)
870  icr1 = ncrnr+(icrn-1)
871  icr1 = mod1(icr1,ncrnr)
872  ivt1 = icface(icr1,ifac)
873  DO 400 ie=1,nelt
874  DO 300 idim=1,ldim
875  side(idim,ifac,ie)=side(idim,ifac,ie)+xyz(idim,ivtx,ie)
876  side( 4,ifac,ie)=side( 4,ifac,ie)+
877  $ ( xyz(idim,ivtx,ie)-xyz(idim,ivt1,ie) )**2
878  300 CONTINUE
879  side(4,ifac,ie)=sqrt( side(4,ifac,ie) )
880  400 CONTINUE
881  500 CONTINUE
882  avwght=1.0/float(ncrnr)
883  CALL cmult(side,avwght,24*nelt)
884 C
885 c call exitt
886  RETURN
887  END
888 c-----------------------------------------------------------------------
889  subroutine verrhe
890 C
891 C 8 Mar 1989 21:58:26 PFF
892 C Verify right-handedness of given elements.
893 C
894  include 'SIZE'
895  include 'INPUT'
896  include 'PARALLEL'
897  include 'SCRCT'
898  include 'TOPOL'
899  LOGICAL IFYES,IFCSTT
900 C
901  ifcstt=.true.
902  IF (.NOT.if3d) THEN
903  DO 1000 ie=1,nelt
904 C
905 C CRSS2D(A,B,O) = (A-O) X (B-O)
906 C
907  c1=crss2d(xyz(1,2,ie),xyz(1,3,ie),xyz(1,1,ie))
908  c2=crss2d(xyz(1,4,ie),xyz(1,1,ie),xyz(1,2,ie))
909  c3=crss2d(xyz(1,1,ie),xyz(1,4,ie),xyz(1,3,ie))
910  c4=crss2d(xyz(1,3,ie),xyz(1,2,ie),xyz(1,4,ie))
911 C
912  IF (c1.LE.0.0.OR.c2.LE.0.0.OR.
913  $ c3.LE.0.0.OR.c4.LE.0.0 ) THEN
914 C
915  ieg=lglel(ie)
916  WRITE(6,800) ieg,c1,c2,c3,c4
917  call exitt
918  800 FORMAT(/,2x,'WARNINGa: Detected non-right-handed element.',
919  $ /,2x,'Number',i8,' C1-4:',4e12.4)
920  ifcstt=.false.
921 C CALL QUERY(IFYES,'Proceed ')
922 C IF (.NOT.IFYES) GOTO 9000
923  ENDIF
924  1000 CONTINUE
925 C
926 C Else 3-D:
927 C
928  ELSE
929  DO 2000 ie=1,nelt
930 C
931 C VOLUM0(A,B,C,O) = (A-O)X(B-O).(C-O)
932 C
933  v1= volum0(xyz(1,2,ie),xyz(1,3,ie),xyz(1,5,ie),xyz(1,1,ie))
934  v2= volum0(xyz(1,4,ie),xyz(1,1,ie),xyz(1,6,ie),xyz(1,2,ie))
935  v3= volum0(xyz(1,1,ie),xyz(1,4,ie),xyz(1,7,ie),xyz(1,3,ie))
936  v4= volum0(xyz(1,3,ie),xyz(1,2,ie),xyz(1,8,ie),xyz(1,4,ie))
937  v5=-volum0(xyz(1,6,ie),xyz(1,7,ie),xyz(1,1,ie),xyz(1,5,ie))
938  v6=-volum0(xyz(1,8,ie),xyz(1,5,ie),xyz(1,2,ie),xyz(1,6,ie))
939  v7=-volum0(xyz(1,5,ie),xyz(1,8,ie),xyz(1,3,ie),xyz(1,7,ie))
940  v8=-volum0(xyz(1,7,ie),xyz(1,6,ie),xyz(1,4,ie),xyz(1,8,ie))
941 C
942  IF (v1.LE.0.0.OR.v2.LE.0.0.OR.
943  $ v3.LE.0.0.OR.v4.LE.0.0.OR.
944  $ v5.LE.0.0.OR.v6.LE.0.0.OR.
945  $ v7.LE.0.0.OR.v8.LE.0.0 ) THEN
946 C
947  ieg=lglel(ie)
948  WRITE(6,1800) ieg,v1,v2,v3,v4,v5,v6,v7,v8
949  call exitt
950  1800 FORMAT(/,2x,'WARNINGb: Detected non-right-handed element.',
951  $ /,2x,'Number',i8,' V1-8:',4e12.4
952  $ /,2x,' ',4x,' ',4e12.4)
953  ifcstt=.false.
954  ENDIF
955  2000 CONTINUE
956  ENDIF
957 C
958  9000 CONTINUE
959 C
960 C Print out results from right-handed check
961 C
962  IF (.NOT.ifcstt) WRITE(6,2001)
963 C
964 C Check consistency accross all processors.
965 C
966  CALL gllog(ifcstt,.false.)
967 C
968  IF (.NOT.ifcstt) THEN
969  IF (nid.EQ.0) WRITE(6,2003) nelgt
970  call exitt
971  ELSE
972  IF (nio.EQ.0) WRITE(6,2002) nelgt
973  ENDIF
974 C
975  2001 FORMAT(//,' Elemental geometry not right-handed, ABORTING'
976  $ ,' in routine VERRHE.')
977  2002 FORMAT(' Right-handed check complete for',i12,' elements. OK.')
978  2003 FORMAT(' Right-handed check failed for',i12,' elements.'
979  $ ,' Exiting in routine VERRHE.')
980  RETURN
981  END
982 c-----------------------------------------------------------------------
983  FUNCTION volum0(P1,P2,P3,P0)
984 C
985 C 3
986 C Given four points in R , (P1,P2,P3,P0), VOLUM0 returns
987 C the volume enclosed by the parallelagram defined by the
988 C vectors { (P1-P0),(P2-P0),(P3-P0) }. This routine has
989 C the nice feature that if the 3 vectors so defined are
990 C not right-handed then the volume returned is negative.
991 C
992  REAL p1(3),p2(3),p3(3),p0(3)
993 C
994  u1=p1(1)-p0(1)
995  u2=p1(2)-p0(2)
996  u3=p1(3)-p0(3)
997 C
998  v1=p2(1)-p0(1)
999  v2=p2(2)-p0(2)
1000  v3=p2(3)-p0(3)
1001 C
1002  w1=p3(1)-p0(1)
1003  w2=p3(2)-p0(2)
1004  w3=p3(3)-p0(3)
1005 C
1006  cross1 = u2*v3-u3*v2
1007  cross2 = u3*v1-u1*v3
1008  cross3 = u1*v2-u2*v1
1009 C
1010  volum0 = w1*cross1 + w2*cross2 + w3*cross3
1011 
1012  RETURN
1013  END
1014 c-----------------------------------------------------------------------
1015  FUNCTION crss2d(XY1,XY2,XY0)
1016  REAL xy1(2),xy2(2),xy0(2)
1017 C
1018  v1x=xy1(1)-xy0(1)
1019  v2x=xy2(1)-xy0(1)
1020  v1y=xy1(2)-xy0(2)
1021  v2y=xy2(2)-xy0(2)
1022  crss2d = v1x*v2y - v1y*v2x
1023 C
1024  RETURN
1025  END
1026 c-----------------------------------------------------------------------
1027  subroutine facind (kx1,kx2,ky1,ky2,kz1,kz2,nx,ny,nz,iface)
1028 c ifcase in preprocessor notation
1029  kx1=1
1030  ky1=1
1031  kz1=1
1032  kx2=nx
1033  ky2=ny
1034  kz2=nz
1035  IF (iface.EQ.1) ky2=1
1036  IF (iface.EQ.2) kx1=nx
1037  IF (iface.EQ.3) ky1=ny
1038  IF (iface.EQ.4) kx2=1
1039  IF (iface.EQ.5) kz2=1
1040  IF (iface.EQ.6) kz1=nz
1041  return
1042  end
1043 c-----------------------------------------------------------------------
1044  subroutine facindr (kx1,kx2,ky1,ky2,kz1,kz2,nx,ny,nz,iface)
1045 
1046 c restricted index set, iface in preprocessor notation
1047 
1048  kx1=2
1049  ky1=2
1050  kz1=2
1051  kx2=nx-1
1052  ky2=ny-1
1053  kz2=nz-1
1054 
1055  if (iface.eq.1) ky1=1
1056  if (iface.eq.1) ky2=1
1057 
1058  if (iface.eq.2) kx1=nx
1059  if (iface.eq.2) kx2=nx
1060 
1061  if (iface.eq.3) ky1=ny
1062  if (iface.eq.3) ky2=ny
1063 
1064  if (iface.eq.4) kx1=1
1065  if (iface.eq.4) kx2=1
1066 
1067  if (iface.eq.5) kz1=1
1068  if (iface.eq.5) kz2=1
1069 
1070  if (iface.eq.6) kz1=nz
1071  if (iface.eq.6) kz2=nz
1072 
1073  return
1074  end
1075 c-----------------------------------------------------------------------
1076  subroutine facev(a,ie,iface,val,nx,ny,nz)
1077 C
1078 C Assign the value VAL to face(IFACE,IE) of array A.
1079 C IFACE is the input in the pre-processor ordering scheme.
1080 C
1081  include 'SIZE'
1082  dimension a(nx,ny,nz,lelt)
1083  CALL facind (kx1,kx2,ky1,ky2,kz1,kz2,nx,ny,nz,iface)
1084  DO 100 iz=kz1,kz2
1085  DO 100 iy=ky1,ky2
1086  DO 100 ix=kx1,kx2
1087  a(ix,iy,iz,ie)=val
1088  100 CONTINUE
1089  RETURN
1090  END
1091 c-----------------------------------------------------------------------
1092  subroutine ifacev(a,ie,iface,val,nx,ny,nz)
1093 C
1094 C Assign the value VAL to face(IFACE,IE) of array A.
1095 C IFACE is the input in the pre-processor ordering scheme.
1096 C
1097  include 'SIZE'
1098  integer a(nx,ny,nz,lelt),val
1099  call facind (kx1,kx2,ky1,ky2,kz1,kz2,nx,ny,nz,iface)
1100  do 100 iz=kz1,kz2
1101  do 100 iy=ky1,ky2
1102  do 100 ix=kx1,kx2
1103  a(ix,iy,iz,ie)=val
1104  100 continue
1105  return
1106  end
1107 c-----------------------------------------------------------------------
1108  subroutine facec(a,b,ie,iface,nx,ny,nz,nel)
1109 C
1110 C Copy the face (IFACE) of B to A.
1111 C IFACE is the input in the pre-processor ordering scheme.
1112 C
1113  dimension a(nx,ny,nz,nel)
1114  dimension b(nx,ny,nz,nel)
1115  CALL facind (kx1,kx2,ky1,ky2,kz1,kz2,nx,ny,nz,iface)
1116  DO 100 iz=kz1,kz2
1117  DO 100 iy=ky1,ky2
1118  DO 100 ix=kx1,kx2
1119  a(ix,iy,iz,ie)=b(ix,iy,iz,ie)
1120  100 CONTINUE
1121  RETURN
1122  END
1123 c-----------------------------------------------------------------------
1124  subroutine combin2(glnm1,glnm2,nglob)
1125 c
1126  write(6,*) 'Hey, who called combin2??? ABORT'
1127  call exitt
1128 c
1129  return
1130  end
1131 c-----------------------------------------------------------------------
1132  subroutine outfldio (x,txt10)
1133  include 'SIZE'
1134  integer x(lx1,ly1,lz1,lelt)
1135  character*10 txt10
1136 C
1137  do ie=1,nelv,2
1138  do iz=1,lz1,1
1139  if (iz.eq.1) write(6,106) txt10,iz,ie
1140  if (iz.gt.1) write(6,107)
1141  i1 = ie+1
1142  do j=ly1,1,-1
1143  write(6,105) (x(i,j,iz,ie),i=1,lx1)
1144  $ , (x(i,j,iz,i1),i=1,lx1)
1145  enddo
1146  enddo
1147  enddo
1148 C
1149  107 FORMAT(' ')
1150  105 FORMAT(4i6,20x,4i6)
1151  106 FORMAT( /,5x,' ^ ',/,
1152  $ 5x,' Y | ',/,
1153  $ 5x,' | ',a10,/,
1154  $ 5x,' +----> ','Plane = ',i2,'/',i2,/,
1155  $ 5x,' X ')
1156 C
1157  return
1158  end
1159 c-----------------------------------------------------------------------
1160  subroutine outfldi (x,txt10)
1161  include 'SIZE'
1162  integer x(lx1,ly1,lz1,lelt)
1163  character*10 txt10
1164 c
1165  character*6 s(20,20)
1166 c
1167  if (lx1.ne.4 .or. nelv.gt.3) return
1168 c
1169  write(6,106) txt10,ie,ie
1170  106 FORMAT( /,5x,' ^ ',/,
1171  $ 5x,' Y | ',/,
1172  $ 5x,' | ',a10,/,
1173  $ 5x,' +----> ','elem. = ',i2,'/',i2,/,
1174  $ 5x,' X ')
1175 C
1176 C
1177  call blank(s,6*20*20)
1178  do ie=1,3
1179  if (ie.eq.1) then
1180  jstart = 1
1181  istart = 1
1182  istride = 3
1183  else
1184  jstart = 2 + lx1
1185  istart = 1
1186  istride = 1
1187  if (ie.eq.2) istart = 7
1188  endif
1189 c
1190  i=istart
1191  do iy=ly1,1,-1
1192  j=jstart
1193  do ix=1,lx1
1194  write(s(i,j),6) x(ix,iy,1,ie)
1195  j=j+1
1196  enddo
1197  i=i+istride
1198  enddo
1199  6 format(20i6)
1200  enddo
1201 c
1202  do i=1,10
1203  write(6,7) (s(i,l),l=1,j-1)
1204  enddo
1205  7 format(20a6)
1206 c
1207  write(6,*)
1208  return
1209  end
1210 c-----------------------------------------------------------------------
1211  subroutine outfldr (x,txt10)
1212  include 'SIZE'
1213  real x(lx1,ly1,lz1,lelt)
1214  character*10 txt10
1215 c
1216  character*6 s(20,20)
1217 c
1218  if (lx1.ne.4 .or. nelv.gt.3) return
1219  write(6,106) txt10,ie,ie
1220  106 FORMAT( /,5x,' ^ ',/,
1221  $ 5x,' Y | ',/,
1222  $ 5x,' | ',a10,/,
1223  $ 5x,' +----> ','elem. = ',i2,'/',i2,/,
1224  $ 5x,' X ')
1225  call blank(s,6*20*20)
1226 C
1227 C
1228  do ie=1,3
1229  if (ie.eq.1) then
1230  jstart = 1
1231  istart = 1
1232  istride = 3
1233  else
1234  jstart = 2 + lx1
1235  istart = 1
1236  istride = 1
1237  if (ie.eq.2) istart = 7
1238  endif
1239 c
1240  i=istart
1241  do iy=ly1,1,-1
1242  j=jstart
1243  do ix=1,lx1
1244  write(s(i,j),6) x(ix,iy,1,ie)
1245  j=j+1
1246  enddo
1247  i=i+istride
1248  enddo
1249  6 format(f6.2)
1250  enddo
1251 c
1252  do i=1,10
1253  write(6,7) (s(i,l),l=1,j-1)
1254  enddo
1255  7 format(20a6)
1256 c
1257  write(6,*)
1258 c
1259  return
1260  end
1261 c-----------------------------------------------------------------------
1262  subroutine checkit(idum)
1263  write(6,*) 'continue?'
1264  read (5,*) idum
1265  return
1266  end
1267 c-----------------------------------------------------------------------
1268  subroutine outfldro (x,txt10,ichk)
1269  include 'SIZE'
1270  include 'TSTEP'
1271  real x(lx1,ly1,lz1,lelt)
1272  character*10 txt10
1273 c
1274  integer idum,e
1275  save idum
1276  data idum /3/
1277  if (idum.lt.0) return
1278 c
1279 C
1280  mtot = lx1*ly1*lz1*nelv
1281  if (lx1.gt.8.or.nelv.gt.16) return
1282  xmin = glmin(x,mtot)
1283  xmax = glmax(x,mtot)
1284 
1285  do ie=1,1
1286  ne = 1
1287  do k=1,1
1288  write(6,116) txt10,k,ie,xmin,xmax,istep,time
1289  write(6,117)
1290  do j=ly1,1,-1
1291  if (lx1.eq.2) write(6,102) ((x(i,j,k,e+1),i=1,lx1),e=0,ne)
1292  if (lx1.eq.3) write(6,103) ((x(i,j,k,e+1),i=1,lx1),e=0,ne)
1293  if (lx1.eq.4) write(6,104) ((x(i,j,k,e+1),i=1,lx1),e=0,ne)
1294  if (lx1.eq.5) write(6,105) ((x(i,j,k,e+1),i=1,lx1),e=0,ne)
1295  if (lx1.eq.6) write(6,106) ((x(i,j,k,e+1),i=1,lx1),e=0,ne)
1296  if (lx1.eq.7) write(6,107) ((x(i,j,k,e+1),i=1,lx1),e=0,ne)
1297  if (lx1.eq.8) write(6,118) ((x(i,j,k,e+1),i=1,lx1),e=0,ne)
1298  enddo
1299  enddo
1300  enddo
1301 
1302  102 FORMAT(4(2f9.5,2x))
1303  103 FORMAT(4(3f9.5,2x))
1304  104 FORMAT(4(4f7.3,2x))
1305  105 FORMAT(5f9.5,10x,5f9.5)
1306  106 FORMAT(6f9.5,5x,6f9.5)
1307  107 FORMAT(7f8.4,5x,7f8.4)
1308  108 FORMAT(8f8.4,4x,8f8.4)
1309  118 FORMAT(8f12.9)
1310 c
1311  116 FORMAT( /,5x,' ^ ',/,
1312  $ 5x,' Y | ',/,
1313  $ 5x,' | ',a10,/,
1314  $ 5x,' +----> ','Plane = ',i2,'/',i2,2x,2e12.4,/,
1315  $ 5x,' X ','Step =',i9,f15.5)
1316  117 FORMAT(' ')
1317 c
1318  if (ichk.eq.1.and.idum.gt.0) call checkit(idum)
1319  return
1320  end
1321 c-----------------------------------------------------------------------
1322  subroutine outfldrv (x,txt10,ichk) ! writes to unit=40+ifield
1323  include 'SIZE'
1324  include 'TSTEP'
1325  real x(lx1,ly1,lz1,lelt)
1326  character*10 txt10
1327 c
1328  integer idum,e
1329  save idum
1330  data idum /3/
1331  if (idum.lt.0) return
1332  m = 40 + ifield
1333 c
1334 C
1335  mtot = lx1*ly1*lz1*nelv
1336  if (lx1.gt.7.or.nelv.gt.16) return
1337  xmin = glmin(x,mtot)
1338  xmax = glmax(x,mtot)
1339 c
1340  rnel = nelv
1341  snel = sqrt(rnel)+.1
1342  ne = snel
1343  ne1 = nelv-ne+1
1344  do ie=ne1,1,-ne
1345  l=ie-1
1346  do k=1,1
1347  if (ie.eq.ne1) write(m,116) txt10,k,ie,xmin,xmax,istep,time
1348  write(m,117)
1349  do j=ly1,1,-1
1350  if (lx1.eq.2) write(m,102) ((x(i,j,k,e+l),i=1,lx1),e=1,ne)
1351  if (lx1.eq.3) write(m,103) ((x(i,j,k,e+l),i=1,lx1),e=1,ne)
1352  if (lx1.eq.4) write(m,104) ((x(i,j,k,e+l),i=1,lx1),e=1,ne)
1353  if (lx1.eq.5) write(m,105) ((x(i,j,k,e+l),i=1,lx1),e=1,ne)
1354  if (lx1.eq.6) write(m,106) ((x(i,j,k,e+l),i=1,lx1),e=1,ne)
1355  if (lx1.eq.7) write(m,107) ((x(i,j,k,e+l),i=1,lx1),e=1,ne)
1356  if (lx1.eq.8) write(m,108) ((x(i,j,k,e+l),i=1,lx1),e=1,ne)
1357  enddo
1358  enddo
1359  enddo
1360 
1361 C
1362  102 FORMAT(4(2f9.5,2x))
1363  103 FORMAT(4(3f9.5,2x))
1364  104 FORMAT(4(4f7.3,2x))
1365  105 FORMAT(5f9.5,10x,5f9.5)
1366  106 FORMAT(6f9.5,5x,6f9.5)
1367  107 FORMAT(7f8.4,5x,7f8.4)
1368  108 FORMAT(8f8.4,4x,8f8.4)
1369 c
1370  116 FORMAT( /,5x,' ^ ',/,
1371  $ 5x,' Y | ',/,
1372  $ 5x,' | ',a10,/,
1373  $ 5x,' +----> ','Plane = ',i2,'/',i2,2x,2e12.4,/,
1374  $ 5x,' X ','Step =',i9,f15.5)
1375  117 FORMAT(' ')
1376 c
1377  if (ichk.eq.1.and.idum.gt.0) call checkit(idum)
1378  return
1379  end
1380 c-----------------------------------------------------------------------
1381  subroutine outfldrv0 (x,txt10,ichk)
1382  include 'SIZE'
1383  include 'TSTEP'
1384  real x(lx1,ly1,lz1,lelt)
1385  character*10 txt10
1386 c
1387  integer idum,e
1388  save idum
1389  data idum /3/
1390  if (idum.lt.0) return
1391 c
1392 C
1393  mtot = lx1*ly1*lz1*nelv
1394  if (lx1.gt.7.or.nelv.gt.16) return
1395  xmin = glmin(x,mtot)
1396  xmax = glmax(x,mtot)
1397 c
1398  rnel = nelv
1399  snel = sqrt(rnel)+.1
1400  ne = snel
1401  ne1 = nelv-ne+1
1402  do ie=ne1,1,-ne
1403  l=ie-1
1404  do k=1,1
1405  if (ie.eq.ne1) write(6,116) txt10,k,ie,xmin,xmax,istep,time
1406  write(6,117)
1407  do j=ly1,1,-1
1408  if (lx1.eq.2) write(6,102) ((x(i,j,k,e+l),i=1,lx1),e=1,ne)
1409  if (lx1.eq.3) write(6,103) ((x(i,j,k,e+l),i=1,lx1),e=1,ne)
1410  if (lx1.eq.4) write(6,104) ((x(i,j,k,e+l),i=1,lx1),e=1,ne)
1411  if (lx1.eq.5) write(6,105) ((x(i,j,k,e+l),i=1,lx1),e=1,ne)
1412  if (lx1.eq.6) write(6,106) ((x(i,j,k,e+l),i=1,lx1),e=1,ne)
1413  if (lx1.eq.7) write(6,107) ((x(i,j,k,e+l),i=1,lx1),e=1,ne)
1414  if (lx1.eq.8) write(6,108) ((x(i,j,k,e+l),i=1,lx1),e=1,ne)
1415  enddo
1416  enddo
1417  enddo
1418 
1419 C
1420  102 FORMAT(4(2f9.5,2x))
1421  103 FORMAT(4(3f9.5,2x))
1422  104 FORMAT(4(4f7.3,2x))
1423  105 FORMAT(5f9.5,10x,5f9.5)
1424  106 FORMAT(6f9.5,5x,6f9.5)
1425  107 FORMAT(7f8.4,5x,7f8.4)
1426  108 FORMAT(8f8.4,4x,8f8.4)
1427 c
1428  116 FORMAT( /,5x,' ^ ',/,
1429  $ 5x,' Y | ',/,
1430  $ 5x,' | ',a10,/,
1431  $ 5x,' +----> ','Plane = ',i2,'/',i2,2x,2e12.4,/,
1432  $ 5x,' X ','Step =',i9,f15.5)
1433  117 FORMAT(' ')
1434 c
1435  if (ichk.eq.1.and.idum.gt.0) call checkit(idum)
1436  return
1437  end
1438 c-----------------------------------------------------------------------
1439  subroutine outfldrp0 (x,txt10,ichk)
1440  include 'SIZE'
1441  include 'TSTEP'
1442  real x(lx2,ly2,lz2,lelt)
1443  character*10 txt10
1444 c
1445  integer idum,e
1446  save idum
1447  data idum /3/
1448  if (idum.lt.0) return
1449 c
1450 C
1451  mtot = lx2*ly2*lz2*nelv
1452  if (lx2.gt.7.or.nelv.gt.16) return
1453  xmin = glmin(x,mtot)
1454  xmax = glmax(x,mtot)
1455 c
1456  rnel = nelv
1457  snel = sqrt(rnel)+.1
1458  ne = snel
1459  ne1 = nelv-ne+1
1460  do ie=ne1,1,-ne
1461  l=ie-1
1462  do k=1,1
1463  if (ie.eq.ne1) write(6,116) txt10,k,ie,xmin,xmax,istep,time
1464  write(6,117)
1465  do j=ly2,1,-1
1466  if (lx2.eq.2) write(6,102) ((x(i,j,k,e+l),i=1,lx2),e=1,ne)
1467  if (lx2.eq.3) write(6,103) ((x(i,j,k,e+l),i=1,lx2),e=1,ne)
1468  if (lx2.eq.4) write(6,104) ((x(i,j,k,e+l),i=1,lx2),e=1,ne)
1469  if (lx2.eq.5) write(6,105) ((x(i,j,k,e+l),i=1,lx2),e=1,ne)
1470  if (lx2.eq.6) write(6,106) ((x(i,j,k,e+l),i=1,lx2),e=1,ne)
1471  if (lx2.eq.7) write(6,107) ((x(i,j,k,e+l),i=1,lx2),e=1,ne)
1472  if (lx2.eq.8) write(6,108) ((x(i,j,k,e+l),i=1,lx2),e=1,ne)
1473  enddo
1474  enddo
1475  enddo
1476 
1477 C
1478  102 FORMAT(4(2f9.5,2x))
1479  103 FORMAT(4(3f9.5,2x))
1480  104 FORMAT(4(4f7.3,2x))
1481  105 FORMAT(5f9.5,10x,5f9.5)
1482  106 FORMAT(6f9.5,5x,6f9.5)
1483  107 FORMAT(7f8.4,5x,7f8.4)
1484  108 FORMAT(8f8.4,4x,8f8.4)
1485 c
1486  116 FORMAT( /,5x,' ^ ',/,
1487  $ 5x,' Y | ',/,
1488  $ 5x,' | ',a10,/,
1489  $ 5x,' +----> ','Plane = ',i2,'/',i2,2x,2e12.4,/,
1490  $ 5x,' X ','Step =',i9,f15.5)
1491  117 FORMAT(' ')
1492 c
1493  if (ichk.eq.1.and.idum.gt.0) call checkit(idum)
1494  return
1495  end
1496 c-----------------------------------------------------------------------
1497  subroutine outfldrp (x,txt10,ichk) ! writes out into unit = 40+ifield
1498  include 'SIZE'
1499  include 'TSTEP'
1500  real x(lx2,ly2,lz2,lelt)
1501  character*10 txt10
1502 c
1503  integer idum,e
1504  save idum
1505  data idum /3/
1506  if (idum.lt.0) return
1507  m = 40 + ifield ! unit #
1508 c
1509 c
1510 C
1511  mtot = lx2*ly2*lz2*nelv
1512  if (lx2.gt.7.or.nelv.gt.16) return
1513  xmin = glmin(x,mtot)
1514  xmax = glmax(x,mtot)
1515 c
1516  rnel = nelv
1517  snel = sqrt(rnel)+.1
1518  ne = snel
1519  ne1 = nelv-ne+1
1520  do ie=ne1,1,-ne
1521  l=ie-1
1522  do k=1,1
1523  if (ie.eq.ne1) write(m,116) txt10,k,ie,xmin,xmax,istep,time
1524  write(6,117)
1525  do j=ly2,1,-1
1526  if (lx2.eq.2) write(m,102) ((x(i,j,k,e+l),i=1,lx2),e=1,ne)
1527  if (lx2.eq.3) write(m,103) ((x(i,j,k,e+l),i=1,lx2),e=1,ne)
1528  if (lx2.eq.4) write(m,104) ((x(i,j,k,e+l),i=1,lx2),e=1,ne)
1529  if (lx2.eq.5) write(m,105) ((x(i,j,k,e+l),i=1,lx2),e=1,ne)
1530  if (lx2.eq.6) write(m,106) ((x(i,j,k,e+l),i=1,lx2),e=1,ne)
1531  if (lx2.eq.7) write(m,107) ((x(i,j,k,e+l),i=1,lx2),e=1,ne)
1532  if (lx2.eq.8) write(m,108) ((x(i,j,k,e+l),i=1,lx2),e=1,ne)
1533  enddo
1534  enddo
1535  enddo
1536 
1537 C
1538  102 FORMAT(4(2f9.5,2x))
1539  103 FORMAT(4(3f9.5,2x))
1540  104 FORMAT(4(4f7.3,2x))
1541  105 FORMAT(5f9.5,10x,5f9.5)
1542  106 FORMAT(6f9.5,5x,6f9.5)
1543  107 FORMAT(7f8.4,5x,7f8.4)
1544  108 FORMAT(8f8.4,4x,8f8.4)
1545 c
1546  116 FORMAT( /,5x,' ^ ',/,
1547  $ 5x,' Y | ',/,
1548  $ 5x,' | ',a10,/,
1549  $ 5x,' +----> ','Plane = ',i2,'/',i2,2x,2e12.4,/,
1550  $ 5x,' X ','Step =',i9,f15.5)
1551  117 FORMAT(' ')
1552 c
1553  if (ichk.eq.1.and.idum.gt.0) call checkit(idum)
1554  return
1555  end
1556 c-----------------------------------------------------------------------
1557  subroutine outmatp(a,m,n,name6,ie)
1558  include 'SIZE'
1559  include 'PARALLEL'
1560  real a(m,n)
1561  character*6 name6
1562 c
1563  if (nelv.gt.1) return
1564  do jid=0,np-1
1565  imax = iglmax(jid,1)
1566  if (jid.eq.nid) then
1567  write(6,*)
1568  write(6,*) nid,ie,' matrix: ',name6,m,n
1569  n12 = min(n,12)
1570  do i=1,m
1571  write(6,6) ie,name6,(a(i,j),j=1,n12)
1572  enddo
1573  6 format(i3,1x,a6,12f9.5)
1574  write(6,*)
1575  call flush_io()
1576  endif
1577  imax = iglmax(jid,1)
1578  enddo
1579  return
1580  end
1581 c-----------------------------------------------------------------------
1582  subroutine gs_chkr(glo_num)
1583  include 'SIZE'
1584  include 'TOTAL'
1585 
1586  integer glo_num(lx1,ly1,lz1,lelt)
1587  integer e,eg
1588 
1589  do ipass = 1,2
1590  iquick = 0
1591  if (nid.eq.0) iquick=glo_num(ipass,1,1,1)
1592  iquick = iglmax(iquick,1)
1593 
1594  do e=1,nelv
1595  do k=1,lz1
1596  do j=1,ly1
1597  do i=1,lx1
1598  if (glo_num(i,j,k,e).eq.iquick) then
1599  eg = lglel(e)
1600  write(6,1) nid,i,j,k,e,eg,iquick,ipass
1601  $ ,xm1(i,j,k,e),ym1(i,j,k,e),zm1(i,j,k,e)
1602  1 format(i12,3i4,2i12,i12,i2,1p3e12.4,' iquick')
1603  endif
1604  enddo
1605  enddo
1606  enddo
1607  enddo
1608  enddo
1609 
1610  return
1611  end
1612 c-----------------------------------------------------------------------
1613  subroutine setup_mesh_dssum ! Set up dssum for mesh
1614 
1615  include 'SIZE'
1616  include 'TOTAL'
1617  include 'NONCON'
1618 
1619  common /c_is1/ glo_num(1*lx1*ly1*lz1*lelv)
1620  integer*8 glo_num
1621  common /ivrtx/ vertex((2**ldim)*lelt)
1622  integer vertex
1623 
1624  parameter(lxyz=lx1*ly1*lz1)
1625  common /scrns/ enum(lxyz,lelt)
1626  $ , rnx(lxyz,lelt) , rny(lxyz,lelt) , rnz(lxyz,lelt)
1627  $ , tnx(lxyz,lelt) , tny(lxyz,lelt) , tnz(lxyz,lelt)
1628  common /scruz/ snx(lxz) , sny(lxz) , snz(lxz) , efc(lxz)
1629  common /scrsf/ jvrtex((2**ldim),lelt)
1630 
1631  integer e,f,eg
1632 
1633 
1634  gsh_fld(0)=gsh_fld(1)
1635  if (iftmsh(0)) gsh_fld(0)=gsh_fld(2)
1636 
1637  ifield = 0
1638  nel = nelfld(0)
1639  nxyz = lx1*ly1*lz1
1640  nxz = lx1*lz1
1641  n = nel*nxyz
1642  nface = 2*ldim
1643 
1644 
1645  iflag=0
1646  do e=1,nel
1647  do f=1,nface
1648  if (cbc(f,e,1).eq.'msi'.or.cbc(f,e,1).eq.'MSI') iflag=1
1649  enddo
1650  enddo
1651  iflag = iglmax(iflag,1)
1652  if (iflag.eq.0) return
1653 
1654 
1655 c We need to differentiate elements according to unit normal on msi face.
1656 
1657 c NOTE that this code assumes we do not have two adjacent faces on a given
1658 c element that both have cbc = msi!
1659 
1660  call rzero(rnx,n)
1661  call rzero(rny,n)
1662  call rzero(rnz,n)
1663  call rzero(tnx,n)
1664  call rzero(tny,n)
1665  call rzero(tnz,n)
1666 
1667  do e=1,nel
1668  re = lglel(e)
1669  call cfill(enum(1,e),re,nxyz)
1670  enddo
1671 
1672  call dsop(enum,'min')
1673 
1674 
1675  do e=1,nel
1676  eg = lglel(e)
1677  do f=1,nface
1678  if (cbc(f,e,1).eq.'msi'.or.cbc(f,e,1).eq.'MSI') then
1679  call facexs (efc,enum(1,e),f,0) ! enum-->efc
1680  do i=1,nxz
1681  jg = efc(i)+0.1
1682  if (jg.eq.eg) then ! this is the controlling face
1683  snx(i) = unx(i,1,f,e)
1684  sny(i) = uny(i,1,f,e)
1685  snz(i) = unz(i,1,f,e)
1686  endif
1687  enddo
1688  call facexv (snx,sny,snz,rnx(1,e),rny(1,e),rnz(1,e),f,1) ! s-->r
1689  call facexv (unx(1,1,f,e),uny(1,1,f,e),unz(1,1,f,e) ! Control
1690  $ ,tnx(1,e),tny(1,e),tnz(1,e),f,1) ! u-->r
1691  endif
1692  enddo
1693  enddo
1694 
1695  call opdssum(rnx,rny,rnz)
1696 
1697  nv = nel*(2**ldim)
1698  call icopy(jvrtex,vertex,nv) ! Save vertex
1699  mvertx=iglmax(jvrtex,nv)
1700 
1701 
1702 c Now, check to see if normal is aligned with incoming normal
1703 
1704  nsx=lx1-1
1705  nsy=ly1-1
1706  nsz=max(1,lz1-1)
1707 
1708  do e=1,nel
1709  l=0
1710  do k=1,lz1,nsz
1711  do j=1,ly1,nsy
1712  do i=1,lx1,nsx
1713  l=l+1
1714  m=i + lx1*(j-1) + lx1*ly1*(k-1)
1715  dot = tnx(m,e)*rnx(m,e)+tny(m,e)*rny(m,e)+tnz(m,e)*rnz(m,e)
1716  if (dot.lt.-.5) jvrtex(l,e)=jvrtex(l,e)+mvertx
1717  enddo
1718  enddo
1719  enddo
1720  enddo
1721 
1722 
1723  call setupds(gsh_fld(0),lx1,ly1,lz1,nelv,nelgv,jvrtex,glo_num)
1724 
1725  return
1726  end
1727 c-----------------------------------------------------------------------
1728  subroutine outfldnx(x,txt10,nx,ny)
1729  include 'SIZE'
1730  real x(nx,ny,4)
1731  character*10 txt10
1732  integer e,g
1733 
1734  write(6,106) txt10,nel,nx
1735  106 FORMAT( /,5x,' ^ ',/,
1736  $ 5x,' Y | ',/,
1737  $ 5x,' | ',a10,/,
1738  $ 5x,' +----> ','elem. = ',i2,'/',i2,/,
1739  $ 5x,' X ')
1740 
1741  do e=3,1,-2
1742  g=e+1
1743  i=istart
1744  do j=ny,1,-1
1745  if (nx.eq.3) write(6,3) (x(i,j,e),i=1,nx),(x(i,j,g),i=1,nx)
1746  if (nx.eq.4) write(6,4) (x(i,j,e),i=1,nx),(x(i,j,g),i=1,nx)
1747  if (nx.eq.5) write(6,5) (x(i,j,e),i=1,nx),(x(i,j,g),i=1,nx)
1748  if (nx.eq.6) write(6,6) (x(i,j,e),i=1,nx),(x(i,j,g),i=1,nx)
1749  if (nx.eq.7) write(6,7) (x(i,j,e),i=1,nx),(x(i,j,g),i=1,nx)
1750  3 format(3f8.4,3x,3f8.4)
1751  4 format(4f8.4,3x,4f8.4)
1752  5 format(5f8.4,3x,5f8.4)
1753  6 format(6f8.4,3x,6f8.4)
1754  7 format(7f8.4,3x,7f8.4)
1755  enddo
1756  write(6,*)
1757  enddo
1758  write(6,*)
1759 
1760  return
1761  end
1762 c-----------------------------------------------------------------------
subroutine setrzer
Definition: bdry.f:145
subroutine chkaxcb
Definition: bdry.f:252
void exitt()
Definition: comm_mpi.f:604
subroutine outfldrp0(x, txt10, ichk)
Definition: connect1.f:1440
subroutine initds
Definition: connect1.f:175
subroutine setedge
Definition: connect1.f:223
subroutine facev(a, ie, iface, val, nx, ny, nz)
Definition: connect1.f:1077
subroutine outfldrv(x, txt10, ichk)
Definition: connect1.f:1323
subroutine setup_mesh_dssum
Definition: connect1.f:1614
subroutine setup_topo
Definition: connect1.f:3
subroutine genxyzl
Definition: connect1.f:716
subroutine setside
Definition: connect1.f:816
subroutine outmatp(a, m, n, name6, ie)
Definition: connect1.f:1558
subroutine outfldrp(x, txt10, ichk)
Definition: connect1.f:1498
subroutine facec(a, b, ie, iface, nx, ny, nz, nel)
Definition: connect1.f:1109
subroutine combin2(glnm1, glnm2, nglob)
Definition: connect1.f:1125
subroutine checkit(idum)
Definition: connect1.f:1263
function crss2d(XY1, XY2, XY0)
Definition: connect1.f:1016
function volum0(P1, P2, P3, P0)
Definition: connect1.f:984
subroutine facind(kx1, kx2, ky1, ky2, kz1, kz2, nx, ny, nz, iface)
Definition: connect1.f:1028
subroutine outfldnx(x, txt10, nx, ny)
Definition: connect1.f:1729
subroutine dsset(nx, ny, nz)
Definition: connect1.f:553
subroutine outfldio(x, txt10)
Definition: connect1.f:1133
subroutine verrhe
Definition: connect1.f:890
subroutine outfldi(x, txt10)
Definition: connect1.f:1161
subroutine outfldrv0(x, txt10, ichk)
Definition: connect1.f:1382
subroutine facindr(kx1, kx2, ky1, ky2, kz1, kz2, nx, ny, nz, iface)
Definition: connect1.f:1045
subroutine ifacev(a, ie, iface, val, nx, ny, nz)
Definition: connect1.f:1093
subroutine verify
Definition: connect1.f:800
subroutine outfldr(x, txt10)
Definition: connect1.f:1212
subroutine gs_chkr(glo_num)
Definition: connect1.f:1583
subroutine outfldro(x, txt10, ichk)
Definition: connect1.f:1269
subroutine dsop(u, op, nx, ny, nz)
Definition: dssum.f:101
subroutine setupds(gs_handle, nx, ny, nz, nel, melg, vertex, glo_num)
Definition: dssum.f:2
subroutine dssum(u, nx, ny, nz)
Definition: dssum.f:34
real function dot(V1, V2, N)
Definition: genxyz.f:885
subroutine get_vert
Definition: map2.f:126
function glmin(a, n)
Definition: math.f:973
function mod1(i, n)
Definition: math.f:509
subroutine icopy(a, b, n)
Definition: math.f:289
subroutine copy(a, b, n)
Definition: math.f:260
subroutine blank(A, N)
Definition: math.f:19
function iglmax(a, n)
Definition: math.f:913
subroutine izero(a, n)
Definition: math.f:215
subroutine gllog(la, lb)
Definition: math.f:986
subroutine cmult(a, const, n)
Definition: math.f:315
subroutine cfill(a, b, n)
Definition: math.f:244
subroutine rone(a, n)
Definition: math.f:230
subroutine invcol1(a, n)
Definition: math.f:62
subroutine rzero(a, n)
Definition: math.f:208
function glmax(a, n)
Definition: math.f:960
subroutine cbcmesh
Definition: mvmesh.f:3
subroutine opdssum(a, b, c)
Definition: navier1.f:2582
Definition: crs_amg.c:1093
subroutine flush_io
Definition: subs1.f:1561
subroutine setcdof
Definition: subs2.f:1742
subroutine facexv(A1, A2, A3, B1, B2, B3, IFACE1, IOP)
Definition: subs2.f:170
subroutine facexs(A, B, IFACE1, IOP)
Definition: subs2.f:125