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