source: CIVL/examples/fortran/nek5000/core/coef.f

main
Last change on this file was ea777aa, checked in by Alex Wilton <awilton@…>, 3 years ago

Moved examples, include, build_default.properties, common.xml, and README out from dev.civl.com into the root of the repo.

git-svn-id: svn://vsl.cis.udel.edu/civl/trunk@5704 fb995dde-84ed-4084-dfe6-e5aef3e2452c

  • Property mode set to 100755
File size: 55.4 KB
Line 
1 subroutine genwz
2C-----------------------------------------------------------------
3C
4C GENERATE
5C
6C - DERIVATIVE OPERATORS
7C - INTERPOLATION OPERATORS
8C - WEIGHTS
9C - COLLOCATION POINTS
10C
11C ASSOCIATED WITH THE
12C
13C - GAUSS-LOBATTO LEGENDRE MESH (SUFFIX M1/M2/M3)
14C - GAUSS LEGENDRE MESH (SUFFIX M2)
15C - GAUSS-LOBATTO JACOBI MESH (SUFFIX M1/M2/M3)
16C
17C-----------------------------------------------------------------
18C
19 INCLUDE 'SIZE'
20 INCLUDE 'WZ'
21 INCLUDE 'DXYZ'
22 INCLUDE 'IXYZ'
23 INCLUDE 'INPUT'
24
25 REAL TMP(LY1,LY1),TMPT(LY1,LY1)
26C
27 IF (ldim.EQ.2) THEN
28C
29C*** Two-dimensional case **********************
30C
31C
32C Gauss-Lobatto Legendre mesh (suffix M1)
33C Generate collocation points and weights
34C
35 CALL ZWGLL (ZGM1(1,1),WXM1,lx1)
36 CALL ZWGLL (ZGM1(1,2),WYM1,ly1)
37 ZGM1(lz1,3) = 0.
38 WZM1(lz1) = 1.
39 DO 100 IY=1,ly1
40 DO 100 IX=1,lx1
41 W3M1(IX,IY,1)=WXM1(IX)*WYM1(IY)
42 100 CONTINUE
43C
44C Compute derivative matrices
45C
46 CALL DGLL (DXM1,DXTM1,ZGM1(1,1),lx1,lx1)
47 CALL DGLL (DYM1,DYTM1,ZGM1(1,2),ly1,ly1)
48 CALL RZERO (DZM1 ,lz1*lz1)
49 CALL RZERO (DZTM1,lz1*lz1)
50C
51C Gauss Legendre mesh (suffix M2)
52C Generate collocation points and weights
53C
54 IF(IFSPLIT)THEN
55 CALL ZWGLL (ZGM2(1,1),WXM2,lx2)
56 CALL ZWGLL (ZGM2(1,2),WYM2,ly2)
57 ELSE
58 CALL ZWGL (ZGM2(1,1),WXM2,lx2)
59 CALL ZWGL (ZGM2(1,2),WYM2,ly2)
60 ENDIF
61 ZGM2(lz2,3) = 0.
62 WZM2(lz2) = 1.
63 DO 200 IY=1,ly2
64 DO 200 IX=1,lx2
65 W3M2(IX,IY,1)=WXM2(IX)*WYM2(IY)
66 200 CONTINUE
67C
68C Gauss-Lobatto Legendre mesh (suffix M3).
69C Generate collocation points and weights.
70C
71 CALL ZWGLL (ZGM3(1,1),WXM3,lx3)
72 CALL ZWGLL (ZGM3(1,2),WYM3,ly3)
73 ZGM3(lz3,3) = 0.
74 WZM3(lz3) = 1.
75 DO 300 IY=1,ly3
76 DO 300 IX=1,lx3
77 W3M3(IX,IY,1)=WXM3(IX)*WYM3(IY)
78 300 CONTINUE
79C
80C Compute derivative matrices
81C
82 CALL DGLL (DXM3,DXTM3,ZGM3(1,1),lx3,lx3)
83 CALL DGLL (DYM3,DYTM3,ZGM3(1,2),ly3,ly3)
84 CALL RZERO (DZM3 ,lz3*lz3)
85 CALL RZERO (DZTM3,lz3*lz3)
86C
87C Generate interpolation operators for the staggered mesh
88C
89 CALL IGLLM (IXM12,IXTM12,ZGM1(1,1),ZGM2(1,1),lx1,lx2,lx1,lx2)
90 CALL IGLLM (IYM12,IYTM12,ZGM1(1,2),ZGM2(1,2),ly1,ly2,ly1,ly2)
91 IZM12 (lz2,lz1) = 1.
92 IZTM12(lz1,lz2) = 1.
93C
94C NOTE: The splitting scheme has only one mesh!!!!!
95C
96 IF (IFSPLIT) THEN
97 CALL IGLLM (IXM21,IXTM21,ZGM1(1,1),ZGM2(1,1),lx1,lx2,lx1,lx2)
98 CALL IGLLM (IYM21,IYTM21,ZGM1(1,2),ZGM2(1,2),ly1,ly2,ly1,ly2)
99 ELSE
100 CALL IGLM (IXM21,IXTM21,ZGM2(1,1),ZGM1(1,1),lx2,lx1,lx2,lx1)
101 CALL IGLM (IYM21,IYTM21,ZGM2(1,2),ZGM1(1,2),ly2,ly1,ly2,ly1)
102 ENDIF
103 IZM21 (lz1,lz2) = 1.
104 IZTM21(lz2,lz1) = 1.
105C
106C Compute derivative operators for the staggered mesh
107C
108 IF(IFSPLIT)THEN
109 CALL COPY (DXM12, DXM1, lx1*lx2)
110 CALL COPY (DXTM12,DXTM1,lx1*lx2)
111 CALL COPY (DYM12, DYM1, ly1*ly2)
112 CALL COPY (DYTM12,DYTM1,ly1*ly2)
113 CALL COPY (DZM12, DZM1, lz1*lz2)
114 CALL COPY (DZTM12,DZTM1,lz1*lz2)
115 ELSE
116 CALL DGLLGL (DXM12,DXTM12,ZGM1(1,1),ZGM2(1,1),IXM12,
117 $ lx1,lx2,lx1,lx2)
118 CALL DGLLGL (DYM12,DYTM12,ZGM1(1,2),ZGM2(1,2),IYM12,
119 $ ly1,ly2,ly1,ly2)
120 DZM12 (lz2,lz1) = 0.
121 DZTM12(lz2,lz1) = 0.
122 ENDIF
123C
124C Compute interpolation operators for the geometry mesh M3.
125C
126 CALL IGLLM (IXM13,IXTM13,ZGM1(1,1),ZGM3(1,1),lx1,lx3,lx1,lx3)
127 CALL IGLLM (IYM13,IYTM13,ZGM1(1,2),ZGM3(1,2),ly1,ly3,ly1,ly3)
128 CALL IGLLM (IXM31,IXTM31,ZGM3(1,1),ZGM1(1,1),lx3,lx1,lx3,lx1)
129 CALL IGLLM (IYM31,IYTM31,ZGM3(1,2),ZGM1(1,2),ly3,ly1,ly3,ly1)
130 IZM13 (lz3,lz1) = 1.
131 IZTM13(lz1,lz3) = 1.
132 IZM31 (lz1,lz3) = 1.
133 IZTM31(lz3,lz1) = 1.
134C
135C
136 IF (IFAXIS) THEN
137C
138C Special treatment for the axisymmetric case
139C Generate additional points, weights, derivative operators and
140C interpolation operators required for elements close to the axis.
141C
142C
143C Gauss-Lobatto Jacobi mesh (suffix M1).
144C Generate collocation points and weights (alpha=0, beta=1).
145C
146 ALPHA = 0.
147 BETA = 1.
148 CALL ZWGLJ (ZAM1,WAM1,ly1,ALPHA,BETA)
149 DO 400 IY=1,ly1
150 DO 400 IX=1,lx1
151 W2AM1(IX,IY)=WXM1(IX)*WAM1(IY)
152 W2CM1(IX,IY)=WXM1(IX)*WYM1(IY)
153 400 CONTINUE
154C
155C Compute derivative matrices
156C
157 CALL COPY (DCM1,DYM1,ly1*ly1)
158 CALL COPY (DCTM1,DYTM1,ly1*ly1)
159 CALL DGLJ (DAM1,DATM1,ZAM1,ly1,ly1,ALPHA,BETA)
160C
161C Gauss Jacobi mesh (suffix M2)
162C Generate collocation points and weights
163C
164 IF(IFSPLIT)THEN
165 CALL ZWGLJ (ZAM2,WAM2,ly2,ALPHA,BETA)
166 ELSE
167 CALL ZWGJ (ZAM2,WAM2,ly2,ALPHA,BETA)
168 ENDIF
169 DO 500 IY=1,ly2
170 DO 500 IX=1,lx2
171 W2CM2(IX,IY)=WXM2(IX)*WYM2(IY)
172 W2AM2(IX,IY)=WXM2(IX)*WAM2(IY)
173 500 CONTINUE
174C
175C Gauss-Lobatto Jacobi mesh (suffix M3).
176C Generate collocation points and weights.
177C
178 CALL ZWGLJ (ZAM3,WAM3,ly3,ALPHA,BETA)
179 DO 600 IY=1,ly3
180 DO 600 IX=1,lx3
181 W2CM3(IX,IY)=WXM3(IX)*WYM3(IY)
182 W2AM3(IX,IY)=WXM3(IX)*WAM3(IY)
183 600 CONTINUE
184C
185C Compute derivative matrices
186C
187 CALL COPY (DCM3,DYM3,ly3*ly3)
188 CALL COPY (DCTM3,DYTM3,ly3*ly3)
189 CALL DGLJ (DAM3,DATM3,ZAM3,ly3,ly3,ALPHA,BETA)
190C
191C Generate interpolation operators for the staggered mesh
192C
193 CALL COPY (ICM12,IYM12,ly2*ly1)
194 CALL COPY (ICTM12,IYTM12,ly1*ly2)
195 CALL IGLJM (IAM12,IATM12,ZAM1,ZAM2,ly1,ly2,ly1,ly2,ALPHA,BETA)
196 CALL COPY (ICM21,IYM21,ly1*ly2)
197 CALL COPY (ICTM21,IYTM21,ly2*ly1)
198 IF (IFSPLIT) THEN
199 CALL IGLJM (IAM21,IATM21,ZAM2,ZAM1,ly1,ly2,ly1,ly2,ALPHA,BETA)
200 ELSE
201 CALL IGJM (IAM21,IATM21,ZAM2,ZAM1,ly2,ly1,ly2,ly1,ALPHA,BETA)
202 ENDIF
203C
204C Compute derivative operators for the staggered mesh
205C
206 CALL COPY (DCM12,DYM12,ly2*ly1)
207 CALL COPY (DCTM12,DYTM12,ly1*ly2)
208 IF(IFSPLIT)THEN
209 CALL COPY (DAM12, DAM1, ly1*ly2)
210 CALL COPY (DATM12,DATM1,ly1*ly2)
211 ELSE
212 CALL DGLJGJ (DAM12,DATM12,ZAM1,ZAM2,IAM12,
213 $ ly1,ly2,ly1,ly2,ALPHA,BETA)
214 ENDIF
215C
216C Compute interpolation operators for the geometry mesh M3.
217C
218 CALL COPY (ICM13,IYM13,ly3*ly1)
219 CALL COPY (ICTM13,IYTM13,ly1*ly3)
220 CALL IGLJM (IAM13,IATM13,ZAM1,ZAM3,ly1,ly3,ly1,ly3,ALPHA,BETA)
221 CALL COPY (ICM31,IYM31,ly1*ly3)
222 CALL COPY (ICTM31,IYTM31,ly3*ly1)
223 CALL IGLJM (IAM31,IATM31,ZAM3,ZAM1,ly3,ly1,ly3,ly1,ALPHA,BETA)
224C
225C Compute interpolation operators between Gauss-Lobatto Jacobi
226C and Gauss-Lobatto Legendre (to be used in PREPOST).
227C
228 CALL IGLJM(IAJL1,IATJL1,ZAM1,ZGM1(1,2),ly1,ly1,ly1,ly1,ALPHA,BETA)
229 IF (IFSPLIT) THEN
230 CALL IGLJM(IAJL2,IATJL2,ZAM2,ZGM2(1,2),ly2,ly2,ly2,ly2,ALPHA,BETA)
231 ELSE
232 CALL IGJM (IAJL2,IATJL2,ZAM2,ZGM2(1,2),ly2,ly2,ly2,ly2,ALPHA,BETA)
233 ENDIF
234
235 CALL INVMT(IAJL1 ,IALJ1 ,TMP ,ly1)
236 CALL INVMT(IATJL1,IATLJ1,TMPT,ly1)
237 CALL MXM (IATJL1,ly1,IATLJ1,ly1,TMPT,ly1)
238 CALL MXM (IAJL1 ,ly1,IALJ1 ,ly1,TMP ,ly1)
239
240C
241C Compute interpolation operators between Gauss-Lobatto Legendre
242C and Gauss-Lobatto Jacobi (to be used in subr. genxyz IN postpre).
243C
244c
245c This call is not right, and these arrays are not used. 3/27/02. pff
246c CALL IGLLM(IALJ3,IATLJ3,ZGM3(1,2),ZAM3,ly3,ly3,ly3,ly3,ALPHA,BETA)
247 CALL IGLJM(IALJ3,IATLJ3,ZGM3(1,2),ZAM3,ly3,ly3,ly3,ly3,ALPHA,BETA)
248C
249 ENDIF
250C
251C
252 ELSE
253C
254C*** Three-dimensional case ************************************
255C
256C
257C Gauss-Lobatto Legendre mesh (suffix M1)
258C Generate collocation points and weights
259C
260 CALL ZWGLL (ZGM1(1,1),WXM1,lx1)
261 CALL ZWGLL (ZGM1(1,2),WYM1,ly1)
262 CALL ZWGLL (ZGM1(1,3),WZM1,lz1)
263 DO 700 IZ=1,lz1
264 DO 700 IY=1,ly1
265 DO 700 IX=1,lx1
266 W3M1(IX,IY,IZ)=WXM1(IX)*WYM1(IY)*WZM1(IZ)
267 700 CONTINUE
268C
269C Compute derivative matrices
270C
271 CALL DGLL (DXM1,DXTM1,ZGM1(1,1),lx1,lx1)
272 CALL DGLL (DYM1,DYTM1,ZGM1(1,2),ly1,ly1)
273 CALL DGLL (DZM1,DZTM1,ZGM1(1,3),lz1,lz1)
274C
275C Gauss Legendre mesh (suffix M2)
276C Generate collocation points and weights
277C
278 IF(IFSPLIT)THEN
279 CALL ZWGLL (ZGM2(1,1),WXM2,lx2)
280 CALL ZWGLL (ZGM2(1,2),WYM2,ly2)
281 CALL ZWGLL (ZGM2(1,3),WZM2,lz2)
282 ELSE
283 CALL ZWGL (ZGM2(1,1),WXM2,lx2)
284 CALL ZWGL (ZGM2(1,2),WYM2,ly2)
285 CALL ZWGL (ZGM2(1,3),WZM2,lz2)
286 ENDIF
287 DO 800 IZ=1,lz2
288 DO 800 IY=1,ly2
289 DO 800 IX=1,lx2
290 W3M2(IX,IY,IZ)=WXM2(IX)*WYM2(IY)*WZM2(IZ)
291 800 CONTINUE
292C
293C Gauss-Loabtto Legendre mesh (suffix M3).
294C Generate collocation points and weights.
295C
296 CALL ZWGLL (ZGM3(1,1),WXM3,lx3)
297 CALL ZWGLL (ZGM3(1,2),WYM3,ly3)
298 CALL ZWGLL (ZGM3(1,3),WZM3,lz3)
299 DO 900 IZ=1,lz3
300 DO 900 IY=1,ly3
301 DO 900 IX=1,lx3
302 W3M3(IX,IY,IZ)=WXM3(IX)*WYM3(IY)*WZM3(IZ)
303 900 CONTINUE
304C
305C Compute derivative matrices
306C
307 CALL DGLL (DXM3,DXTM3,ZGM3(1,1),lx3,lx3)
308 CALL DGLL (DYM3,DYTM3,ZGM3(1,2),ly3,ly3)
309 CALL DGLL (DZM3,DZTM3,ZGM3(1,3),lz3,lz3)
310C
311C Generate interpolation operators for the staggered mesh
312C
313 CALL IGLLM (IXM12,IXTM12,ZGM1(1,1),ZGM2(1,1),lx1,lx2,lx1,lx2)
314 CALL IGLLM (IYM12,IYTM12,ZGM1(1,2),ZGM2(1,2),ly1,ly2,ly1,ly2)
315 CALL IGLLM (IZM12,IZTM12,ZGM1(1,3),ZGM2(1,3),lz1,lz2,lz1,lz2)
316C
317C NOTE: The splitting scheme has only one mesh!!!!!
318C
319 IF (IFSPLIT) THEN
320 CALL IGLLM (IXM21,IXTM21,ZGM1(1,1),ZGM2(1,1),lx1,lx2,lx1,lx2)
321 CALL IGLLM (IYM21,IYTM21,ZGM1(1,2),ZGM2(1,2),ly1,ly2,ly1,ly2)
322 CALL IGLLM (IZM21,IZTM21,ZGM1(1,3),ZGM2(1,3),lz1,lz2,lz1,lz2)
323 ELSE
324 CALL IGLM (IXM21,IXTM21,ZGM2(1,1),ZGM1(1,1),lx2,lx1,lx2,lx1)
325 CALL IGLM (IYM21,IYTM21,ZGM2(1,2),ZGM1(1,2),ly2,ly1,ly2,ly1)
326 CALL IGLM (IZM21,IZTM21,ZGM2(1,3),ZGM1(1,3),lz2,lz1,lz2,lz1)
327 ENDIF
328C
329C Compute derivative operators for the staggered mesh
330C
331 IF(IFSPLIT)THEN
332 CALL COPY (DXM12, DXM1, lx1*lx2)
333 CALL COPY (DXTM12,DXTM1,lx1*lx2)
334 CALL COPY (DYM12, DYM1, ly1*ly2)
335 CALL COPY (DYTM12,DYTM1,ly1*ly2)
336 CALL COPY (DZM12, DZM1, lz1*lz2)
337 CALL COPY (DZTM12,DZTM1,lz1*lz2)
338 ELSE
339 CALL DGLLGL (DXM12,DXTM12,ZGM1(1,1),ZGM2(1,1),IXM12,
340 $ lx1,lx2,lx1,lx2)
341 CALL DGLLGL (DYM12,DYTM12,ZGM1(1,2),ZGM2(1,2),IYM12,
342 $ ly1,ly2,ly1,ly2)
343 CALL DGLLGL (DZM12,DZTM12,ZGM1(1,3),ZGM2(1,3),IZM12,
344 $ lz1,lz2,lz1,lz2)
345 ENDIF
346C
347C Compute interpolation operators for the geometry mesh M3.
348C
349 CALL IGLLM (IXM13,IXTM13,ZGM1(1,1),ZGM3(1,1),lx1,lx3,lx1,lx3)
350 CALL IGLLM (IYM13,IYTM13,ZGM1(1,2),ZGM3(1,2),ly1,ly3,ly1,ly3)
351 CALL IGLLM (IZM13,IZTM13,ZGM1(1,3),ZGM3(1,3),lz1,lz3,lz1,lz3)
352 CALL IGLLM (IXM31,IXTM31,ZGM3(1,1),ZGM1(1,1),lx3,lx1,lx3,lx1)
353 CALL IGLLM (IYM31,IYTM31,ZGM3(1,2),ZGM1(1,2),ly3,ly1,ly3,ly1)
354 CALL IGLLM (IZM31,IZTM31,ZGM3(1,3),ZGM1(1,3),lz3,lz1,lz3,lz1)
355C
356 ENDIF
357C
358 RETURN
359 END
360 subroutine geom1 (xm3,ym3,zm3)
361C-----------------------------------------------------------------------
362C
363C Routine to generate all elemental geometric data for mesh 1.
364C
365C Velocity formulation : global-to-local mapping based on mesh 3
366C Stress formulation : global-to-local mapping based on mesh 1
367C
368C-----------------------------------------------------------------------
369 INCLUDE 'SIZE'
370 INCLUDE 'GEOM'
371 INCLUDE 'INPUT'
372 INCLUDE 'TSTEP'
373C
374C Note : XM3,YM3,ZM3 should come from COMMON /SCRUZ/.
375C
376 DIMENSION XM3(LX3,LY3,LZ3,1)
377 $ , YM3(LX3,LY3,LZ3,1)
378 $ , ZM3(LX3,LY3,LZ3,1)
379C
380 IF (IFGMSH3 .AND. ISTEP.EQ.0) THEN
381 CALL GLMAPM3 (XM3,YM3,ZM3)
382 ELSE
383 CALL GLMAPM1
384 ENDIF
385C
386 CALL GEODAT1
387C
388 RETURN
389 END
390 subroutine glmapm3 (xm3,ym3,zm3)
391C-------------------------------------------------------------------
392C
393C Routine to generate mapping data based on mesh 3
394C (Gauss-Legendre Lobatto meshes).
395C
396C XRM3, YRM3, ZRM3 - dx/dr, dy/dr, dz/dr
397C XSM3, YSM3, ZSM3 - dx/ds, dy/ds, dz/ds
398C XTM3, YTM3, ZTM3 - dx/dt, dy/dt, dz/dt
399C RXM3, RYM3, RZM3 - dr/dx, dr/dy, dr/dz
400C SXM3, SYM3, SZM3 - ds/dx, ds/dy, ds/dz
401C TXM3, TYM3, TZM3 - dt/dx, dt/dy, dt/dz
402C JACM3 - Jacobian
403C
404C------------------------------------------------------------------
405 INCLUDE 'SIZE'
406 INCLUDE 'TOTAL'
407C
408C Note : work arrays for mesh 3 in scratch commons will be
409C changed after exit of routine.
410C
411 COMMON /SCRNS/ XRM3 (LX3,LY3,LZ3,LELT)
412 $ , XSM3 (LX3,LY3,LZ3,LELT)
413 $ , XTM3 (LX3,LY3,LZ3,LELT)
414 $ , YRM3 (LX3,LY3,LZ3,LELT)
415 $ , YSM3 (LX3,LY3,LZ3,LELT)
416 $ , YTM3 (LX3,LY3,LZ3,LELT)
417 $ , ZRM3 (LX3,LY3,LZ3,LELT)
418 COMMON /CTMP0/ ZSM3 (LX3,LY3,LZ3,LELT)
419 $ , ZTM3 (LX3,LY3,LZ3,LELT)
420 COMMON /CTMP1/ RXM3 (LX3,LY3,LZ3,LELT)
421 $ , RYM3 (LX3,LY3,LZ3,LELT)
422 $ , RZM3 (LX3,LY3,LZ3,LELT)
423 $ , SXM3 (LX3,LY3,LZ3,LELT)
424 COMMON /SCRMG/ SYM3 (LX3,LY3,LZ3,LELT)
425 $ , SZM3 (LX3,LY3,LZ3,LELT)
426 $ , TXM3 (LX3,LY3,LZ3,LELT)
427 $ , TYM3 (LX3,LY3,LZ3,LELT)
428 COMMON /SCREV/ TZM3 (LX3,LY3,LZ3,LELT)
429 $ , JACM3(LX3,LY3,LZ3,LELT)
430 REAL JACM3
431 DIMENSION XM3(LX3,LY3,LZ3,1)
432 $ , YM3(LX3,LY3,LZ3,1)
433 $ , ZM3(LX3,LY3,LZ3,1)
434C
435C
436 NXY3 = lx3*ly3
437 NYZ3 = ly3*lz3
438 NXYZ3 = lx3*ly3*lz3
439 NTOT3 = NXYZ3*NELT
440 NXYZ1 = lx1*ly1*lz1
441 NTOT1 = NXYZ1*NELT
442C
443C
444C Compute isoparametric partials.
445C
446
447 IF (ldim.EQ.2) THEN
448C
449C Two-dimensional case
450C
451 DO 200 IEL=1,NELT
452C
453C Use the appropriate derivative- and interpolation operator in
454C the y-direction (= radial direction if axisymmetric).
455C
456 IF (IFAXIS) THEN
457 ly33 = ly3*ly3
458 IF (IFRZER(IEL)) THEN
459 CALL COPY (DYTM3,DATM3,ly33)
460 ELSE
461 CALL COPY (DYTM3,DCTM3,ly33)
462 ENDIF
463 ENDIF
464C
465 CALL MXM(DXM3,lx3,XM3(1,1,1,IEL),lx3,XRM3(1,1,1,IEL),ly3)
466 CALL MXM(DXM3,lx3,YM3(1,1,1,IEL),lx3,YRM3(1,1,1,IEL),ly3)
467 CALL MXM(XM3(1,1,1,IEL),lx3,DYTM3,ly3,XSM3(1,1,1,IEL),ly3)
468 CALL MXM(YM3(1,1,1,IEL),lx3,DYTM3,ly3,YSM3(1,1,1,IEL),ly3)
469C
470 200 CONTINUE
471C
472 CALL RZERO (JACM3,NTOT3)
473 CALL ADDCOL3 (JACM3,XRM3,YSM3,NTOT3)
474 CALL SUBCOL3 (JACM3,XSM3,YRM3,NTOT3)
475C
476 CALL COPY (RXM3,YSM3,NTOT3)
477 CALL COPY (RYM3,XSM3,NTOT3)
478 CALL CHSIGN (RYM3,NTOT3)
479 CALL COPY (SXM3,YRM3,NTOT3)
480 CALL CHSIGN (SXM3,NTOT3)
481 CALL COPY (SYM3,XRM3,NTOT3)
482C
483 ELSE
484C
485C Three-dimensional case
486C
487 DO 300 IEL=1,NELT
488C
489 CALL MXM(DXM3,lx3,XM3(1,1,1,IEL),lx3,XRM3(1,1,1,IEL),NYZ3)
490 CALL MXM(DXM3,lx3,YM3(1,1,1,IEL),lx3,YRM3(1,1,1,IEL),NYZ3)
491 CALL MXM(DXM3,lx3,ZM3(1,1,1,IEL),lx3,ZRM3(1,1,1,IEL),NYZ3)
492C
493 DO 310 IZ=1,lz3
494 CALL MXM(XM3(1,1,IZ,IEL),lx3,DYTM3,ly3,XSM3(1,1,IZ,IEL),ly3)
495 CALL MXM(YM3(1,1,IZ,IEL),lx3,DYTM3,ly3,YSM3(1,1,IZ,IEL),ly3)
496 CALL MXM(ZM3(1,1,IZ,IEL),lx3,DYTM3,ly3,ZSM3(1,1,IZ,IEL),ly3)
497 310 CONTINUE
498C
499 CALL MXM(XM3(1,1,1,IEL),NXY3,DZTM3,lz3,XTM3(1,1,1,IEL),lz3)
500 CALL MXM(YM3(1,1,1,IEL),NXY3,DZTM3,lz3,YTM3(1,1,1,IEL),lz3)
501 CALL MXM(ZM3(1,1,1,IEL),NXY3,DZTM3,lz3,ZTM3(1,1,1,IEL),lz3)
502C
503 300 CONTINUE
504C
505 CALL RZERO (JACM3,NTOT3)
506 CALL ADDCOL4 (JACM3,XRM3,YSM3,ZTM3,NTOT3)
507 CALL ADDCOL4 (JACM3,XTM3,YRM3,ZSM3,NTOT3)
508 CALL ADDCOL4 (JACM3,XSM3,YTM3,ZRM3,NTOT3)
509 CALL SUBCOL4 (JACM3,XRM3,YTM3,ZSM3,NTOT3)
510 CALL SUBCOL4 (JACM3,XSM3,YRM3,ZTM3,NTOT3)
511 CALL SUBCOL4 (JACM3,XTM3,YSM3,ZRM3,NTOT3)
512C
513 CALL ASCOL5 (RXM3,YSM3,ZTM3,YTM3,ZSM3,NTOT3)
514 CALL ASCOL5 (RYM3,XTM3,ZSM3,XSM3,ZTM3,NTOT3)
515 CALL ASCOL5 (RZM3,XSM3,YTM3,XTM3,YSM3,NTOT3)
516 CALL ASCOL5 (SXM3,YTM3,ZRM3,YRM3,ZTM3,NTOT3)
517 CALL ASCOL5 (SYM3,XRM3,ZTM3,XTM3,ZRM3,NTOT3)
518 CALL ASCOL5 (SZM3,XTM3,YRM3,XRM3,YTM3,NTOT3)
519 CALL ASCOL5 (TXM3,YRM3,ZSM3,YSM3,ZRM3,NTOT3)
520 CALL ASCOL5 (TYM3,XSM3,ZRM3,XRM3,ZSM3,NTOT3)
521 CALL ASCOL5 (TZM3,XRM3,YSM3,XSM3,YRM3,NTOT3)
522C
523 ENDIF
524C
525C Mapping from space P(n-2) to space P(n) (mesh M3 to mesh M1).
526C
527 IF (ldim.EQ.2) THEN
528 CALL RZERO (RZM1,NTOT1)
529 CALL RZERO (SZM1,NTOT1)
530 CALL RONE (TZM1,NTOT1)
531 ENDIF
532C
533 kerr = 0
534 DO 400 ie=1,NELT
535
536c write(6,*) 'chkj1'
537c call outxm3j(xm3,ym3,jacm3)
538
539 CALL CHKJAC(JACM3(1,1,1,ie),NXYZ3,ie,xm3(1,1,1,ie),
540 $ ym3(1,1,1,ie),zm3(1,1,1,ie),ldim,ierr)
541 if (ierr.eq.1) kerr = kerr+1
542 CALL MAP31 (RXM1(1,1,1,ie),RXM3(1,1,1,ie),ie)
543 CALL MAP31 (RYM1(1,1,1,ie),RYM3(1,1,1,ie),ie)
544 CALL MAP31 (SXM1(1,1,1,ie),SXM3(1,1,1,ie),ie)
545 CALL MAP31 (SYM1(1,1,1,ie),SYM3(1,1,1,ie),ie)
546 IF (ldim.EQ.3) THEN
547 CALL MAP31 (RZM1(1,1,1,ie),RZM3(1,1,1,ie),ie)
548 CALL MAP31 (SZM1(1,1,1,ie),SZM3(1,1,1,ie),ie)
549 CALL MAP31 (TXM1(1,1,1,ie),TXM3(1,1,1,ie),ie)
550 CALL MAP31 (TYM1(1,1,1,ie),TYM3(1,1,1,ie),ie)
551 CALL MAP31 (TZM1(1,1,1,ie),TZM3(1,1,1,ie),ie)
552 ENDIF
553 CALL MAP31 (JACM1(1,1,1,ie),JACM3(1,1,1,ie),ie)
554 CALL MAP31 (XM1(1,1,1,ie),XM3(1,1,1,ie),ie)
555 CALL MAP31 (YM1(1,1,1,ie),YM3(1,1,1,ie),ie)
556 CALL MAP31 (ZM1(1,1,1,ie),ZM3(1,1,1,ie),ie)
557 400 CONTINUE
558 kerr = iglsum(kerr,1)
559 if (kerr.gt.0) then
560 ifxyo = .true.
561 ifvo = .false.
562 ifpo = .false.
563 ifto = .false.
564 param(66) = 4
565 call outpost(vx,vy,vz,pr,t,'xyz')
566 if (nid.eq.0) write(6,*) 'Jac error 3, setting p66=4, ifxyo=t'
567 call exitt
568 endif
569
570 call invers2(jacmi,jacm1,ntot1)
571
572 RETURN
573 END
574 subroutine glmapm1
575C-----------------------------------------------------------------------
576C
577C Routine to generate mapping data based on mesh 1
578C (Gauss-Legendre Lobatto meshes).
579C
580C XRM1, YRM1, ZRM1 - dx/dr, dy/dr, dz/dr
581C XSM1, YSM1, ZSM1 - dx/ds, dy/ds, dz/ds
582C XTM1, YTM1, ZTM1 - dx/dt, dy/dt, dz/dt
583C RXM1, RYM1, RZM1 - dr/dx, dr/dy, dr/dz
584C SXM1, SYM1, SZM1 - ds/dx, ds/dy, ds/dz
585C TXM1, TYM1, TZM1 - dt/dx, dt/dy, dt/dz
586C JACM1 - Jacobian
587C
588C-----------------------------------------------------------------------
589 INCLUDE 'SIZE'
590 INCLUDE 'GEOM'
591 INCLUDE 'INPUT'
592 INCLUDE 'SOLN'
593C
594C Note: Subroutines GLMAPM1, GEODAT1, AREA2, SETWGTR and AREA3
595C share the same array structure in Scratch Common /SCRNS/.
596C
597 COMMON /SCRNS/ XRM1(LX1,LY1,LZ1,LELT)
598 $ , YRM1(LX1,LY1,LZ1,LELT)
599 $ , XSM1(LX1,LY1,LZ1,LELT)
600 $ , YSM1(LX1,LY1,LZ1,LELT)
601 $ , XTM1(LX1,LY1,LZ1,LELT)
602 $ , YTM1(LX1,LY1,LZ1,LELT)
603 $ , ZRM1(LX1,LY1,LZ1,LELT)
604 COMMON /CTMP1/ ZSM1(LX1,LY1,LZ1,LELT)
605 $ , ZTM1(LX1,LY1,LZ1,LELT)
606C
607 NXY1 = lx1*ly1
608 NYZ1 = ly1*lz1
609 NXYZ1 = lx1*ly1*lz1
610 NTOT1 = NXYZ1*NELT
611C
612 CALL XYZRST (XRM1,YRM1,ZRM1,XSM1,YSM1,ZSM1,XTM1,YTM1,ZTM1,
613 $ IFAXIS)
614C
615 IF (ldim.EQ.2) THEN
616 CALL RZERO (JACM1,NTOT1)
617 CALL ADDCOL3 (JACM1,XRM1,YSM1,NTOT1)
618 CALL SUBCOL3 (JACM1,XSM1,YRM1,NTOT1)
619 CALL COPY (RXM1,YSM1,NTOT1)
620 CALL COPY (RYM1,XSM1,NTOT1)
621 CALL CHSIGN (RYM1,NTOT1)
622 CALL COPY (SXM1,YRM1,NTOT1)
623 CALL CHSIGN (SXM1,NTOT1)
624 CALL COPY (SYM1,XRM1,NTOT1)
625 CALL RZERO (RZM1,NTOT1)
626 CALL RZERO (SZM1,NTOT1)
627 CALL RONE (TZM1,NTOT1)
628 ELSE
629 CALL RZERO (JACM1,NTOT1)
630 CALL ADDCOL4 (JACM1,XRM1,YSM1,ZTM1,NTOT1)
631 CALL ADDCOL4 (JACM1,XTM1,YRM1,ZSM1,NTOT1)
632 CALL ADDCOL4 (JACM1,XSM1,YTM1,ZRM1,NTOT1)
633 CALL SUBCOL4 (JACM1,XRM1,YTM1,ZSM1,NTOT1)
634 CALL SUBCOL4 (JACM1,XSM1,YRM1,ZTM1,NTOT1)
635 CALL SUBCOL4 (JACM1,XTM1,YSM1,ZRM1,NTOT1)
636 CALL ASCOL5 (RXM1,YSM1,ZTM1,YTM1,ZSM1,NTOT1)
637 CALL ASCOL5 (RYM1,XTM1,ZSM1,XSM1,ZTM1,NTOT1)
638 CALL ASCOL5 (RZM1,XSM1,YTM1,XTM1,YSM1,NTOT1)
639 CALL ASCOL5 (SXM1,YTM1,ZRM1,YRM1,ZTM1,NTOT1)
640 CALL ASCOL5 (SYM1,XRM1,ZTM1,XTM1,ZRM1,NTOT1)
641 CALL ASCOL5 (SZM1,XTM1,YRM1,XRM1,YTM1,NTOT1)
642 CALL ASCOL5 (TXM1,YRM1,ZSM1,YSM1,ZRM1,NTOT1)
643 CALL ASCOL5 (TYM1,XSM1,ZRM1,XRM1,ZSM1,NTOT1)
644 CALL ASCOL5 (TZM1,XRM1,YSM1,XSM1,YRM1,NTOT1)
645 ENDIF
646C
647 kerr = 0
648 DO 500 ie=1,NELT
649 CALL CHKJAC(JACM1(1,1,1,ie),NXYZ1,ie,xm1(1,1,1,ie),
650 $ ym1(1,1,1,ie),zm1(1,1,1,ie),ldim,ierr)
651 if (ierr.ne.0) kerr = kerr+1
652 500 CONTINUE
653 kerr = iglsum(kerr,1)
654 if (kerr.gt.0) then
655 ifxyo = .true.
656 ifvo = .false.
657 ifpo = .false.
658 ifto = .false.
659 param(66) = 4
660 call outpost(vx,vy,vz,pr,t,'xyz')
661 if (nid.eq.0) write(6,*) 'Jac error 1, setting p66=4, ifxyo=t'
662 call exitt
663 endif
664
665 call invers2(jacmi,jacm1,ntot1)
666
667 RETURN
668 END
669 subroutine geodat1
670C-----------------------------------------------------------------------
671C
672C Routine to generate elemental geometric matrices on mesh 1
673C (Gauss-Legendre Lobatto mesh).
674C
675C-----------------------------------------------------------------------
676 INCLUDE 'SIZE'
677 INCLUDE 'GEOM'
678 INCLUDE 'INPUT'
679 INCLUDE 'MASS'
680 INCLUDE 'TSTEP'
681 INCLUDE 'WZ'
682C
683C Note: Subroutines GLMAPM1, GEODAT1, AREA2, SETWGTR and AREA3
684C share the same array structure in Scratch Common /SCRNS/.
685C
686 COMMON /SCRNS/ XRM1(LX1,LY1,LZ1,LELT)
687 $ , YRM1(LX1,LY1,LZ1,LELT)
688 $ , XSM1(LX1,LY1,LZ1,LELT)
689 $ , YSM1(LX1,LY1,LZ1,LELT)
690 $ , XTM1(LX1,LY1,LZ1,LELT)
691 $ , YTM1(LX1,LY1,LZ1,LELT)
692 $ , ZRM1(LX1,LY1,LZ1,LELT)
693 COMMON /CTMP1/ ZSM1(LX1,LY1,LZ1,LELT)
694 $ , ZTM1(LX1,LY1,LZ1,LELT)
695 $ , WJ (LX1,LY1,LZ1,LELT)
696C
697 NXYZ1 = lx1*ly1*lz1
698 NTOT1 = NXYZ1*NELT
699C
700 IF (IFGMSH3 .AND. ISTEP.EQ.0)
701 $ CALL XYZRST (XRM1,YRM1,ZRM1,XSM1,YSM1,ZSM1,XTM1,YTM1,ZTM1,
702 $ IFAXIS)
703C
704 IF (.NOT.IFAXIS) THEN
705 CALL INVERS2 (WJ,JACM1,NTOT1)
706 ELSE
707 DO 500 IEL=1,NELT
708 IF (IFRZER(IEL)) THEN
709 DO 510 J=1,ly1
710 DO 510 I=1,lx1
711 IF (J.GT.1) THEN
712 WJ(I,J,1,IEL) = YM1(I,J,1,IEL)/
713 $ (JACM1(I,J,1,IEL)*(1.+ZAM1(J)))
714 ELSE
715 WJ(I,J,1,IEL) = YSM1(I,J,1,IEL)/JACM1(I,J,1,IEL)
716 ENDIF
717 510 CONTINUE
718 ELSE
719 CALL INVCOL3 (WJ(1,1,1,IEL),YM1(1,1,1,IEL),
720 $ JACM1(1,1,1,IEL),NXYZ1)
721 ENDIF
722 500 CONTINUE
723 ENDIF
724C
725C Compute geometric factors for integrated del-squared operator.
726C
727 IF (ldim.EQ.2) THEN
728 CALL VDOT2 (G1M1,RXM1,RYM1,RXM1,RYM1,NTOT1)
729 CALL VDOT2 (G2M1,SXM1,SYM1,SXM1,SYM1,NTOT1)
730 CALL VDOT2 (G4M1,RXM1,RYM1,SXM1,SYM1,NTOT1)
731 CALL COL2 (G1M1,WJ,NTOT1)
732 CALL COL2 (G2M1,WJ,NTOT1)
733 CALL COL2 (G4M1,WJ,NTOT1)
734 CALL RZERO (G3M1,NTOT1)
735 CALL RZERO (G5M1,NTOT1)
736 CALL RZERO (G6M1,NTOT1)
737 ELSE
738 CALL VDOT3 (G1M1,RXM1,RYM1,RZM1,RXM1,RYM1,RZM1,NTOT1)
739 CALL VDOT3 (G2M1,SXM1,SYM1,SZM1,SXM1,SYM1,SZM1,NTOT1)
740 CALL VDOT3 (G3M1,TXM1,TYM1,TZM1,TXM1,TYM1,TZM1,NTOT1)
741 CALL VDOT3 (G4M1,RXM1,RYM1,RZM1,SXM1,SYM1,SZM1,NTOT1)
742 CALL VDOT3 (G5M1,RXM1,RYM1,RZM1,TXM1,TYM1,TZM1,NTOT1)
743 CALL VDOT3 (G6M1,SXM1,SYM1,SZM1,TXM1,TYM1,TZM1,NTOT1)
744 CALL COL2 (G1M1,WJ,NTOT1)
745 CALL COL2 (G2M1,WJ,NTOT1)
746 CALL COL2 (G3M1,WJ,NTOT1)
747 CALL COL2 (G4M1,WJ,NTOT1)
748 CALL COL2 (G5M1,WJ,NTOT1)
749 CALL COL2 (G6M1,WJ,NTOT1)
750 ENDIF
751C
752C Multiply the geometric factors GiM1,i=1,5 with the
753C weights on mesh M1.
754C
755 DO 580 IEL=1,NELT
756 IF (IFAXIS) CALL SETAXW1 ( IFRZER(IEL) )
757 CALL COL2 (G1M1(1,1,1,IEL),W3M1,NXYZ1)
758 CALL COL2 (G2M1(1,1,1,IEL),W3M1,NXYZ1)
759 CALL COL2 (G4M1(1,1,1,IEL),W3M1,NXYZ1)
760 IF (ldim.EQ.3) THEN
761 CALL COL2 (G3M1(1,1,1,IEL),W3M1,NXYZ1)
762 CALL COL2 (G5M1(1,1,1,IEL),W3M1,NXYZ1)
763 CALL COL2 (G6M1(1,1,1,IEL),W3M1,NXYZ1)
764 ENDIF
765 580 CONTINUE
766C
767C Compute the mass matrix on mesh M1.
768C
769 DO 700 IEL=1,NELT
770 IF (IFAXIS) CALL SETAXW1 ( IFRZER(IEL) )
771 CALL COL3 (BM1 (1,1,1,IEL),JACM1(1,1,1,IEL),W3M1,NXYZ1)
772 IF (IFAXIS) THEN
773 CALL COL3(BAXM1(1,1,1,IEL),JACM1(1,1,1,IEL),W3M1,NXYZ1)
774 IF (IFRZER(IEL)) THEN
775 DO 600 J=1,ly1
776 IF (J.GT.1) THEN
777 DO 610 I=1,lx1
778 BM1(I,J,1,IEL) = BM1(I,J,1,IEL)*YM1(I,J,1,IEL)
779 $ /(1.+ZAM1(J))
780 BAXM1(I,J,1,IEL)=BAXM1(I,J,1,IEL)/(1.+ZAM1(J))
781 610 CONTINUE
782 ELSE
783 DO 620 I=1,lx1
784 BM1(I,J,1,IEL) = BM1(I,J,1,IEL)*YSM1(I,J,1,IEL)
785 BAXM1(I,J,1,IEL)=BAXM1(I,J,1,IEL)
786 620 CONTINUE
787 ENDIF
788 600 CONTINUE
789 ELSE
790 CALL COL2 (BM1(1,1,1,IEL),YM1(1,1,1,IEL),NXYZ1)
791 ENDIF
792 ENDIF
793C
794 700 CONTINUE
795
796 IF(IFAXIS) THEN
797 DO IEL=1,NELT
798 IF(IFRZER(IEL)) THEN
799 DO J=1,ly1
800 DO I=1,lx1
801 IF(J.EQ.1) THEN
802 YINVM1(I,J,1,IEL)=1.0D0/YSM1(I,J,1,IEL)
803 ELSE
804 YINVM1(I,J,1,IEL)=1.0D0/YM1 (I,J,1,IEL)
805 ENDIF
806 ENDDO
807 ENDDO
808 ELSE
809 CALL INVERS2(YINVM1(1,1,1,IEL),YM1(1,1,1,IEL),NXYZ1)
810 ENDIF
811 ENDDO
812 ELSE
813 CALL CFILL(YINVM1,1.0D0,NXYZ1*NELT)
814 ENDIF
815C
816C Compute normals, tangents, and areas on elemental surfaces
817C
818 CALL SETAREA
819C
820 RETURN
821 END
822 subroutine geom2
823C-------------------------------------------------------------------
824C
825C Routine to generate all elemental geometric data for mesh 2
826C (Gauss-Legendre mesh).
827C
828C RXM2, RYM2, RZM2 - dr/dx, dr/dy, dr/dz
829C SXM2, SYM2, SZM2 - ds/dx, ds/dy, ds/dz
830C TXM2, TYM2, TZM2 - dt/dx, dt/dy, dt/dz
831C JACM2 - Jacobian
832C BM2 - Mass matrix
833C
834C------------------------------------------------------------------
835 INCLUDE 'SIZE'
836 INCLUDE 'TOTAL'
837C
838 NXYZ2 = lx2*ly2*lz2
839 NTOT2 = NXYZ2*NELV
840C
841 IF (IFSPLIT) THEN
842C
843C Mesh 1 and 2 are identical
844C
845 CALL COPY (RXM2,RXM1,NTOT2)
846 CALL COPY (RYM2,RYM1,NTOT2)
847 CALL COPY (RZM2,RZM1,NTOT2)
848 CALL COPY (SXM2,SXM1,NTOT2)
849 CALL COPY (SYM2,SYM1,NTOT2)
850 CALL COPY (SZM2,SZM1,NTOT2)
851 CALL COPY (TXM2,TXM1,NTOT2)
852 CALL COPY (TYM2,TYM1,NTOT2)
853 CALL COPY (TZM2,TZM1,NTOT2)
854 CALL COPY (JACM2,JACM1,NTOT2)
855 CALL COPY (BM2,BM1,NTOT2)
856
857 CALL COPY (XM2,XM1,NTOT2)
858 CALL COPY (YM2,YM1,NTOT2)
859 CALL COPY (ZM2,ZM1,NTOT2)
860
861 ELSE
862C
863C Consistent approximation spaces (UZAWA)
864C
865 IF (ldim.EQ.2) THEN
866 CALL RZERO (RZM2,NTOT2)
867 CALL RZERO (SZM2,NTOT2)
868 CALL RONE (TZM2,NTOT2)
869 ENDIF
870C
871 DO 1000 IEL=1,NELV
872C
873C Mapping from mesh M1 to mesh M2
874C
875 CALL MAP12 (RXM2(1,1,1,IEL),RXM1(1,1,1,IEL),IEL)
876 CALL MAP12 (RYM2(1,1,1,IEL),RYM1(1,1,1,IEL),IEL)
877 CALL MAP12 (SXM2(1,1,1,IEL),SXM1(1,1,1,IEL),IEL)
878 CALL MAP12 (SYM2(1,1,1,IEL),SYM1(1,1,1,IEL),IEL)
879 IF (ldim.EQ.3) THEN
880 CALL MAP12 (RZM2(1,1,1,IEL),RZM1(1,1,1,IEL),IEL)
881 CALL MAP12 (SZM2(1,1,1,IEL),SZM1(1,1,1,IEL),IEL)
882 CALL MAP12 (TXM2(1,1,1,IEL),TXM1(1,1,1,IEL),IEL)
883 CALL MAP12 (TYM2(1,1,1,IEL),TYM1(1,1,1,IEL),IEL)
884 CALL MAP12 (TZM2(1,1,1,IEL),TZM1(1,1,1,IEL),IEL)
885 ENDIF
886 CALL MAP12 (JACM2(1,1,1,IEL),JACM1(1,1,1,IEL),IEL)
887C
888 CALL MAP12 (XM2(1,1,1,IEL),XM1(1,1,1,IEL),IEL)
889 CALL MAP12 (YM2(1,1,1,IEL),YM1(1,1,1,IEL),IEL)
890 CALL MAP12 (ZM2(1,1,1,IEL),ZM1(1,1,1,IEL),IEL)
891C
892C Compute the mass matrix on mesh M2.
893C
894 IF (IFAXIS) CALL SETAXW2 ( IFRZER(IEL) )
895 CALL COL3 (BM2(1,1,1,IEL),W3M2,JACM2(1,1,1,IEL),NXYZ2)
896C
897 IF (IFAXIS.AND.IFRZER(IEL)) THEN
898 DO 300 J=1,ly2
899 DO 300 I=1,lx2
900 BM2(I,J,1,IEL) = BM2(I,J,1,IEL)*YM2(I,J,1,IEL)
901 $ /(1.+ZAM2(J))
902 300 CONTINUE
903 ELSEIF (IFAXIS.AND.(.NOT.IFRZER(IEL))) THEN
904 CALL COL2 (BM2(1,1,1,IEL),YM2(1,1,1,IEL),NXYZ2)
905 ENDIF
906 1000 CONTINUE
907C
908 ENDIF
909C
910C Compute inverse of mesh 2 mass matrix, pff 3/5/92
911 CALL INVERS2(BM2INV,BM2,NTOT2)
912C
913 RETURN
914 END
915 subroutine xyzrst (xrm1,yrm1,zrm1,xsm1,ysm1,zsm1,
916 $ XTM1,YTM1,ZTM1,IFAXIS)
917C-----------------------------------------------------------------------
918C
919C Compute global-to-local derivatives on mesh 1.
920C
921C-----------------------------------------------------------------------
922 INCLUDE 'SIZE'
923 INCLUDE 'GEOM'
924 INCLUDE 'DXYZ'
925C
926 DIMENSION XRM1(LX1,LY1,LZ1,1),YRM1(LX1,LY1,LZ1,1)
927 $ , ZRM1(LX1,LY1,LZ1,1),XSM1(LX1,LY1,LZ1,1)
928 $ , YSM1(LX1,LY1,LZ1,1),ZSM1(LX1,LY1,LZ1,1)
929 $ , XTM1(LX1,LY1,LZ1,1),YTM1(LX1,LY1,LZ1,1)
930 $ , ZTM1(LX1,LY1,LZ1,1)
931 LOGICAL IFAXIS
932C
933 NXY1=lx1*ly1
934 NYZ1=ly1*lz1
935C
936 DO 100 IEL=1,NELT
937C
938 IF (IFAXIS) CALL SETAXDY ( IFRZER(IEL) )
939C
940 CALL MXM (DXM1,lx1,XM1(1,1,1,IEL),lx1,XRM1(1,1,1,IEL),NYZ1)
941 CALL MXM (DXM1,lx1,YM1(1,1,1,IEL),lx1,YRM1(1,1,1,IEL),NYZ1)
942 CALL MXM (DXM1,lx1,ZM1(1,1,1,IEL),lx1,ZRM1(1,1,1,IEL),NYZ1)
943C
944 DO 10 IZ=1,lz1
945 CALL MXM (XM1(1,1,IZ,IEL),lx1,DYTM1,ly1,XSM1(1,1,IZ,IEL),ly1)
946 CALL MXM (YM1(1,1,IZ,IEL),lx1,DYTM1,ly1,YSM1(1,1,IZ,IEL),ly1)
947 CALL MXM (ZM1(1,1,IZ,IEL),lx1,DYTM1,ly1,ZSM1(1,1,IZ,IEL),ly1)
948 10 CONTINUE
949C
950 IF (ldim.EQ.3) THEN
951 CALL MXM (XM1(1,1,1,IEL),NXY1,DZTM1,lz1,XTM1(1,1,1,IEL),lz1)
952 CALL MXM (YM1(1,1,1,IEL),NXY1,DZTM1,lz1,YTM1(1,1,1,IEL),lz1)
953 CALL MXM (ZM1(1,1,1,IEL),NXY1,DZTM1,lz1,ZTM1(1,1,1,IEL),lz1)
954 ELSE
955 CALL RZERO (XTM1(1,1,1,IEL),NXY1)
956 CALL RZERO (YTM1(1,1,1,IEL),NXY1)
957 CALL RONE (ZTM1(1,1,1,IEL),NXY1)
958 ENDIF
959C
960 100 CONTINUE
961C
962 RETURN
963 END
964 subroutine chkjac(jac,n,iel,X,Y,Z,ND,IERR)
965c
966 include 'SIZE'
967 include 'PARALLEL'
968C
969C Check the array JAC for a change in sign.
970C
971 REAL JAC(N),x(1),y(1),z(1)
972c
973 ierr = 1
974 SIGN = JAC(1)
975 DO 100 I=2,N
976 IF (SIGN*JAC(I).LE.0.0) THEN
977 ieg = lglel(iel)
978 WRITE(6,101) nid,I,ieg
979 write(6,*) jac(i-1),jac(i)
980 if (ldim.eq.3) then
981 write(6,7) nid,x(i-1),y(i-1),z(i-1)
982 write(6,7) nid,x(i),y(i),z(i)
983 else
984 write(6,7) nid,x(i-1),y(i-1)
985 write(6,7) nid,x(i),y(i)
986 endif
987 7 format(i5,' xyz:',1p3e14.5)
988c if (np.eq.1) call out_xyz_el(x,y,z,iel)
989c ierr=0
990 return
991 ENDIF
992 100 CONTINUE
993 101 FORMAT(//,i5,2x
994 $ ,'ERROR: Vanishing Jacobian near',i7,'th node of element'
995 $ ,I10,'.')
996c
997c
998 ierr = 0
999 RETURN
1000 END
1001c-----------------------------------------------------------------------
1002 subroutine volume
1003C
1004C Compute the volume based on mesh M1 and mesh M2
1005C
1006
1007 include 'SIZE'
1008 include 'ESOLV'
1009 include 'INPUT'
1010 include 'MASS'
1011 include 'TSTEP'
1012 integer e
1013C
1014 volvm1=glsum(bm1,lx1*ly1*lz1*nelv)
1015 volvm2=glsum(bm2,lx2*ly2*lz2*nelv)
1016 voltm1=glsum(bm1,lx1*ly1*lz1*nelt)
1017 voltm2=glsum(bm2,lx2*ly2*lz2*nelt)
1018 mfield=1
1019 if (ifmvbd) mfield=0
1020 nfldt = nfield
1021 if (ifmhd) nfldt = nfield+1
1022
1023 do ifld=mfield,nfldt
1024 if (iftmsh(ifld)) then
1025 volfld(ifld) = voltm1
1026 else
1027 volfld(ifld) = volvm1
1028 endif
1029 enddo
1030
1031c if (nio.eq.0) write(6,*) 'vol_t,vol_v:',voltm1,volvm1
1032
1033
1034 nxyz = lx1*ly1*lz1
1035 do e=1,nelt
1036 volel(e) = vlsum(bm1(1,1,1,e),nxyz)
1037 enddo
1038
1039 return
1040 end
1041c-----------------------------------------------------------------------
1042 subroutine setarea
1043C
1044C Compute surface data: areas, normals and tangents
1045C
1046 INCLUDE 'SIZE'
1047 INCLUDE 'GEOM'
1048 INCLUDE 'INPUT'
1049C
1050 NSRF = 6*lx1*lz1*NELT
1051C
1052 CALL RZERO (AREA,NSRF)
1053 CALL RZERO3 (UNX,UNY,UNZ,NSRF)
1054 CALL RZERO3 (T1X,T1Y,T1Z,NSRF)
1055 CALL RZERO3 (T2X,T2Y,T2Z,NSRF)
1056C
1057 IF (ldim.EQ.2) THEN
1058 CALL AREA2
1059 ELSE
1060 CALL AREA3
1061 ENDIF
1062C
1063 RETURN
1064 END
1065 subroutine area2
1066C--------------------------------------------------------------------
1067C
1068C Compute areas, normals and tangents (2D and Axisymmetric geom.)
1069C
1070C--------------------------------------------------------------------
1071 INCLUDE 'SIZE'
1072 INCLUDE 'GEOM'
1073C
1074C Note: Subroutines GLMAPM1, GEODAT1, AREA2, SETWGTR and AREA3
1075C share the same array structure in Scratch Common /SCRNS/.
1076C
1077 COMMON /SCRNS/ XRM1(LX1,LY1,LZ1,LELT)
1078 $ , YRM1(LX1,LY1,LZ1,LELT)
1079 $ , XSM1(LX1,LY1,LZ1,LELT)
1080 $ , YSM1(LX1,LY1,LZ1,LELT)
1081 COMMON /CTMP0/ WGTR1(LX1,LELT)
1082 $ , WGTR2(LY1,LELT)
1083 $ , WGTR3(LX1,LELT)
1084 $ , WGTR4(LY1,LELT)
1085C
1086 CALL SETWGTR (WGTR1,WGTR2,WGTR3,WGTR4)
1087C
1088C "R"
1089C
1090 DO 100 IEL=1,NELT
1091 DO 100 IY=1,ly1
1092 XS2 = XSM1(lx1,IY,1,IEL)
1093 YS2 = YSM1(lx1,IY,1,IEL)
1094 XS4 = XSM1( 1,IY,1,IEL)
1095 YS4 = YSM1( 1,IY,1,IEL)
1096 SS2 = SQRT( XS2**2 + YS2**2 )
1097 SS4 = SQRT( XS4**2 + YS4**2 )
1098 T1X (IY,1,2,IEL) = XS2 / SS2
1099 T1Y (IY,1,2,IEL) = YS2 / SS2
1100 T1X (IY,1,4,IEL) = -XS4 / SS4
1101 T1Y (IY,1,4,IEL) = -YS4 / SS4
1102 UNX (IY,1,2,IEL) = T1Y(IY,1,2,IEL)
1103 UNY (IY,1,2,IEL) = -T1X(IY,1,2,IEL)
1104 UNX (IY,1,4,IEL) = T1Y(IY,1,4,IEL)
1105 UNY (IY,1,4,IEL) = -T1X(IY,1,4,IEL)
1106 AREA(IY,1,2,IEL) = SS2 * WGTR2(IY,IEL)
1107 AREA(IY,1,4,IEL) = SS4 * WGTR4(IY,IEL)
1108 100 CONTINUE
1109C
1110C "S"
1111C
1112 DO 200 IEL=1,NELT
1113 DO 200 IX=1,lx1
1114 XR1 = XRM1(IX, 1,1,IEL)
1115 YR1 = YRM1(IX, 1,1,IEL)
1116 XR3 = XRM1(IX,ly1,1,IEL)
1117 YR3 = YRM1(IX,ly1,1,IEL)
1118 RR1 = SQRT( XR1**2 + YR1**2 )
1119 RR3 = SQRT( XR3**2 + YR3**2 )
1120 T1X (IX,1,1,IEL) = XR1 / RR1
1121 T1Y (IX,1,1,IEL) = YR1 / RR1
1122 T1X (IX,1,3,IEL) = -XR3 / RR3
1123 T1Y (IX,1,3,IEL) = -YR3 / RR3
1124 UNX (IX,1,1,IEL) = T1Y(IX,1,1,IEL)
1125 UNY (IX,1,1,IEL) = -T1X(IX,1,1,IEL)
1126 UNX (IX,1,3,IEL) = T1Y(IX,1,3,IEL)
1127 UNY (IX,1,3,IEL) = -T1X(IX,1,3,IEL)
1128 AREA(IX,1,1,IEL) = RR1 * WGTR1(IX,IEL)
1129 AREA(IX,1,3,IEL) = RR3 * WGTR3(IX,IEL)
1130 200 CONTINUE
1131C
1132 RETURN
1133 END
1134 subroutine setwgtr (wgtr1,wgtr2,wgtr3,wgtr4)
1135C
1136 INCLUDE 'SIZE'
1137 INCLUDE 'GEOM'
1138 INCLUDE 'INPUT'
1139 INCLUDE 'WZ'
1140C
1141C Note: Subroutines GLMAPM1, GEODAT1, AREA2, SETWGTR and AREA3
1142C share the same array structure in Scratch Common /SCRNS/.
1143C
1144 COMMON /SCRNS/ XRM1(LX1,LY1,LZ1,LELT)
1145 $ , YRM1(LX1,LY1,LZ1,LELT)
1146 $ , XSM1(LX1,LY1,LZ1,LELT)
1147 $ , YSM1(LX1,LY1,LZ1,LELT)
1148C
1149 DIMENSION WGTR1(LX1,1)
1150 $ , WGTR2(LY1,1)
1151 $ , WGTR3(LX1,1)
1152 $ , WGTR4(LY1,1)
1153C
1154 IF (IFAXIS) THEN
1155 DO 100 IEL=1,NELT
1156 DO 120 IX=1,lx1
1157 WGTR1(IX,IEL) = YM1(IX, 1,1,IEL) * WXM1(IX)
1158 WGTR3(IX,IEL) = YM1(IX,ly1,1,IEL) * WXM1(IX)
1159 120 CONTINUE
1160 IF ( IFRZER(IEL) ) THEN
1161 IY = 1
1162 WGTR2(IY,IEL) = YSM1(lx1,IY,1,IEL) * WAM1(IY)
1163 WGTR4(IY,IEL) = YSM1( 1,IY,1,IEL) * WAM1(IY)
1164 DO 160 IY=2,ly1
1165 DNR = 1. + ZAM1(IY)
1166 WGTR2(IY,IEL) = YM1(lx1,IY,1,IEL) * WAM1(IY) / DNR
1167 WGTR4(IY,IEL) = YM1( 1,IY,1,IEL) * WAM1(IY) / DNR
1168 160 CONTINUE
1169 ELSE
1170 DO 180 IY=1,ly1
1171 WGTR2(IY,IEL) = YM1(lx1,IY,1,IEL) * WYM1(IY)
1172 WGTR4(IY,IEL) = YM1( 1,IY,1,IEL) * WYM1(IY)
1173 180 CONTINUE
1174 ENDIF
1175 100 CONTINUE
1176 ELSE
1177 DO 200 IEL=1,NELT
1178 CALL COPY (WGTR1(1,IEL),WXM1,lx1)
1179 CALL COPY (WGTR2(1,IEL),WYM1,ly1)
1180 CALL COPY (WGTR3(1,IEL),WXM1,lx1)
1181 CALL COPY (WGTR4(1,IEL),WYM1,ly1)
1182 200 CONTINUE
1183 ENDIF
1184C
1185 RETURN
1186 END
1187 subroutine area3
1188C--------------------------------------------------------------------
1189C
1190C Compute areas, normals and tangents (3D geom.)
1191C
1192C--------------------------------------------------------------------
1193 INCLUDE 'SIZE'
1194 INCLUDE 'WZ'
1195 INCLUDE 'GEOM'
1196C
1197C Note: Subroutines GLMAPM1, GEODAT1, AREA2, SETWGTR and AREA3
1198C share the same array structure in Scratch Common /SCRNS/.
1199C
1200 COMMON /SCRNS/ XRM1(LX1,LY1,LZ1,LELT)
1201 $ , YRM1(LX1,LY1,LZ1,LELT)
1202 $ , XSM1(LX1,LY1,LZ1,LELT)
1203 $ , YSM1(LX1,LY1,LZ1,LELT)
1204 $ , XTM1(LX1,LY1,LZ1,LELT)
1205 $ , YTM1(LX1,LY1,LZ1,LELT)
1206 $ , ZRM1(LX1,LY1,LZ1,LELT)
1207 COMMON /CTMP1/ ZSM1(LX1,LY1,LZ1,LELT)
1208 $ , ZTM1(LX1,LY1,LZ1,LELT)
1209 $ , A (LX1,LY1,LZ1,LELT)
1210 $ , B (LX1,LY1,LZ1,LELT)
1211 COMMON /CTMP0/ C (LX1,LY1,LZ1,LELT)
1212 $ , DOT(LX1,LY1,LZ1,LELT)
1213C
1214 NXY1 = lx1*ly1
1215 NFACE = 2*ldim
1216 NTOT = lx1*ly1*lz1*NELT
1217 NSRF = 6*lx1*ly1*NELT
1218C
1219C "R"
1220C
1221 CALL VCROSS(A,B,C,XSM1,YSM1,ZSM1,XTM1,YTM1,ZTM1,NTOT)
1222 CALL VDOT3 (DOT,A,B,C,A,B,C,NTOT)
1223C
1224 DO 100 IEL=1,NELT
1225 DO 100 IZ=1,lz1
1226 DO 100 IY=1,ly1
1227 WEIGHT = WYM1(IY)*WZM1(IZ)
1228 AREA(IY,IZ,2,IEL) = SQRT(DOT(lx1,IY,IZ,IEL))*WEIGHT
1229 AREA(IY,IZ,4,IEL) = SQRT(DOT( 1,IY,IZ,IEL))*WEIGHT
1230 UNX (IY,IZ,4,IEL) = -A( 1,IY,IZ,IEL)
1231 UNX (IY,IZ,2,IEL) = A(lx1,IY,IZ,IEL)
1232 UNY (IY,IZ,4,IEL) = -B( 1,IY,IZ,IEL)
1233 UNY (IY,IZ,2,IEL) = B(lx1,IY,IZ,IEL)
1234 UNZ (IY,IZ,4,IEL) = -C( 1,IY,IZ,IEL)
1235 UNZ (IY,IZ,2,IEL) = C(lx1,IY,IZ,IEL)
1236 100 CONTINUE
1237C
1238C "S"
1239C
1240 CALL VCROSS(A,B,C,XRM1,YRM1,ZRM1,XTM1,YTM1,ZTM1,NTOT)
1241 CALL VDOT3 (DOT,A,B,C,A,B,C,NTOT)
1242 DO 200 IEL=1,NELT
1243 DO 200 IZ=1,lz1
1244 DO 200 IX=1,lx1
1245 WEIGHT=WXM1(IX)*WZM1(IZ)
1246 AREA(IX,IZ,1,IEL) = SQRT(DOT(IX, 1,IZ,IEL))*WEIGHT
1247 AREA(IX,IZ,3,IEL) = SQRT(DOT(IX,ly1,IZ,IEL))*WEIGHT
1248 UNX (IX,IZ,1,IEL) = A(IX, 1,IZ,IEL)
1249 UNX (IX,IZ,3,IEL) = -A(IX,ly1,IZ,IEL)
1250 UNY (IX,IZ,1,IEL) = B(IX, 1,IZ,IEL)
1251 UNY (IX,IZ,3,IEL) = -B(IX,ly1,IZ,IEL)
1252 UNZ (IX,IZ,1,IEL) = C(IX, 1,IZ,IEL)
1253 UNZ (IX,IZ,3,IEL) = -C(IX,ly1,IZ,IEL)
1254 200 CONTINUE
1255C
1256C "T"
1257C
1258 CALL VCROSS(A,B,C,XRM1,YRM1,ZRM1,XSM1,YSM1,ZSM1,NTOT)
1259 CALL VDOT3 (DOT,A,B,C,A,B,C,NTOT)
1260 DO 300 IEL=1,NELT
1261 DO 300 IX=1,lx1
1262 DO 300 IY=1,ly1
1263 WEIGHT=WXM1(IX)*WYM1(IY)
1264 AREA(IX,IY,5,IEL) = SQRT(DOT(IX,IY, 1,IEL))*WEIGHT
1265 AREA(IX,IY,6,IEL) = SQRT(DOT(IX,IY,lz1,IEL))*WEIGHT
1266 UNX (IX,IY,5,IEL) = -A(IX,IY, 1,IEL)
1267 UNX (IX,IY,6,IEL) = A(IX,IY,lz1,IEL)
1268 UNY (IX,IY,5,IEL) = -B(IX,IY, 1,IEL)
1269 UNY (IX,IY,6,IEL) = B(IX,IY,lz1,IEL)
1270 UNZ (IX,IY,5,IEL) = -C(IX,IY, 1,IEL)
1271 UNZ (IX,IY,6,IEL) = C(IX,IY,lz1,IEL)
1272 300 CONTINUE
1273C
1274 CALL UNITVEC (UNX,UNY,UNZ,NSRF)
1275C
1276C COMPUTE UNIT TANGENT T1
1277C
1278 DO 600 IEL=1,NELT
1279 DO 600 IFC=1,NFACE
1280 IF (IFC.EQ.1 .OR. IFC.EQ.6) THEN
1281 CALL FACEXV (T1X(1,1,IFC,IEL),T1Y(1,1,IFC,IEL),
1282 $ T1Z(1,1,IFC,IEL),
1283 $ XRM1(1,1,1,IEL),YRM1(1,1,1,IEL),
1284 $ ZRM1(1,1,1,IEL),IFC,0)
1285 ELSEIF (IFC.EQ.2 .OR. IFC.EQ.5) THEN
1286 CALL FACEXV (T1X(1,1,IFC,IEL),T1Y(1,1,IFC,IEL),
1287 $ T1Z(1,1,IFC,IEL),
1288 $ XSM1(1,1,1,IEL),YSM1(1,1,1,IEL),
1289 $ ZSM1(1,1,1,IEL),IFC,0)
1290 ELSE
1291 CALL FACEXV (T1X(1,1,IFC,IEL),T1Y(1,1,IFC,IEL),
1292 $ T1Z(1,1,IFC,IEL),
1293 $ XTM1(1,1,1,IEL),YTM1(1,1,1,IEL),
1294 $ ZTM1(1,1,1,IEL),IFC,0)
1295 ENDIF
1296 600 CONTINUE
1297C
1298 CALL UNITVEC (T1X,T1Y,T1Z,NSRF)
1299C
1300C COMPUTE UNIT TANGENT T2 ( T2 = Normal X T1 )
1301C
1302 DO 700 IEL=1,NELT
1303 DO 700 IFC=1,NFACE
1304 CALL VCROSS (T2X(1,1,IFC,IEL),T2Y(1,1,IFC,IEL),
1305 $ T2Z(1,1,IFC,IEL),
1306 $ UNX(1,1,IFC,IEL),UNY(1,1,IFC,IEL),
1307 $ UNZ(1,1,IFC,IEL),
1308 $ T1X(1,1,IFC,IEL),T1Y(1,1,IFC,IEL),
1309 $ T1Z(1,1,IFC,IEL),NXY1)
1310 700 CONTINUE
1311C
1312 RETURN
1313 END
1314 subroutine lagmass
1315C--------------------------------------------------------------------
1316C
1317C Lag the mass matrix (matrices)
1318C
1319C--------------------------------------------------------------------
1320 INCLUDE 'SIZE'
1321 INCLUDE 'MASS'
1322 INCLUDE 'TSTEP'
1323C
1324 NTOT1 = lx1*ly1*lz1*NELT
1325 DO 100 ILAG=NBDINP-1,2,-1
1326 CALL COPY (BM1LAG(1,1,1,1,ILAG),BM1LAG(1,1,1,1,ILAG-1),NTOT1)
1327 100 CONTINUE
1328 CALL COPY (BM1LAG(1,1,1,1,1),BM1,NTOT1)
1329C
1330 RETURN
1331 END
1332C
1333 subroutine setinvm
1334C--------------------------------------------------------------------
1335C
1336C Invert the mass matrix.
1337C
1338C 1) Copy BM1 to BINVM1
1339C 2) Perform direct stiffness summation on BINVM1
1340C 3) Compute BINVM1 = 1/BINVM1
1341C 4) Two inverse mass matrices required because of difference
1342C in DSSUM routine for IMESH=1 and IMESH=2.
1343C
1344C--------------------------------------------------------------------
1345 INCLUDE 'SIZE'
1346 INCLUDE 'MASS'
1347 INCLUDE 'GEOM'
1348 INCLUDE 'INPUT'
1349 INCLUDE 'TSTEP'
1350 INCLUDE 'WZ'
1351
1352 nxyz1 = lx1*ly1*lz1
1353
1354 ifld = ifield
1355
1356csk IF (IFFLOW) THEN ! Velocity mass matrix
1357 IFIELD = 1
1358 NTOT = NXYZ1*NELV
1359 CALL COPY (BINVM1,BM1,NTOT)
1360 CALL DSSUM (BINVM1,lx1,ly1,lz1)
1361 CALL INVCOL1 (BINVM1,NTOT)
1362csk ENDIF
1363
1364
1365 IF (IFHEAT) THEN ! Temperature mass matrix
1366 IFIELD = 2
1367 NTOT = NXYZ1*NELT
1368 CALL COPY (BINTM1,BM1,NTOT)
1369 CALL DSSUM (BINTM1,lx1,ly1,lz1)
1370 CALL INVCOL1 (BINTM1,NTOT)
1371 ENDIF
1372
1373 ifield = ifld
1374
1375 return
1376 end
1377
1378c-----------------------------------------------------------------------
1379
1380 subroutine maprs(y,x,xa,nrest,iel)
1381C
1382C Map the elemental array X from Restart mesh to Y on mesh M1
1383C Conforming elements, i.e. lx1=ly1=lz1.
1384C
1385C---------------------------------------------------------------
1386C
1387 INCLUDE 'SIZE'
1388 INCLUDE 'GEOM'
1389 INCLUDE 'IXYZ'
1390 INCLUDE 'WZ'
1391 INCLUDE 'INPUT'
1392C
1393 REAL X(NREST,NREST,NREST)
1394 REAL Y(LX1,LY1,LZ1)
1395C
1396 REAL XA(lx1,NREST,NREST)
1397 COMMON /CTMP0/ XB(LX1,LY1,LZ1)
1398C
1399 REAL IXRES(LX1,LX1),IXTRES(LX1,LX1)
1400 REAL IYRES(LY1,LY1),IYTRES(LY1,LY1)
1401 REAL IZRES(LZ1,LZ1),IZTRES(LZ1,LZ1)
1402 REAL ZCRES(20),WCRES(20)
1403 REAL ZARES(20),WARES(20)
1404C
1405 NZREST = NREST
1406 IF(lz1.EQ.1) NZREST=1
1407 NYZRES = NREST*NZREST
1408 NXY1 = lx1 *ly1
1409C
1410 CALL ZWGLL (ZCRES,WCRES,NREST)
1411 CALL IGLLM (IXRES,IXTRES,ZCRES,ZGM1,NREST,lx1,NREST,lx1)
1412 IF (.NOT.IFAXIS) THEN
1413 CALL COPY (IYRES,IXRES,lx1*NREST)
1414 CALL COPY (IYTRES,IXTRES,lx1*NREST)
1415 CALL COPY (IZRES,IXRES,lx1*NREST)
1416 CALL COPY (IZTRES,IXTRES,lx1*NREST)
1417 ELSE
1418C
1419C Use the appropriate derivative- and interpolation operator in
1420C the y-direction (= radial direction if axisymmetric).
1421C
1422 IF (IFRZER(IEL)) THEN
1423 ALPHA = 0.
1424 BETA = 1.
1425 CALL ZWGLJ (ZARES,WARES,NREST,ALPHA,BETA)
1426 CALL IGLJM (IYRES,IYTRES,ZARES,ZGM1,NREST,ly1,NREST,ly1,
1427 $ ALPHA,BETA)
1428 ly1R = ly1*NREST
1429 ELSE
1430 CALL COPY (IYRES,IXRES,lx1*NREST)
1431 CALL COPY (IYTRES,IXTRES,lx1*NREST)
1432 ENDIF
1433 ENDIF
1434C
1435 IF (ldim.EQ.2) THEN
1436 CALL MXM (IXRES,lx1,X,NREST,XA,NREST)
1437 CALL MXM (XA,lx1,IYTRES,NREST,Y,ly1)
1438 ELSE
1439 CALL MXM (IXRES,lx1,X,NREST,XA,NYZRES)
1440 DO 100 IZ=1,NZREST
1441 CALL MXM (XA(1,1,IZ),lx1,IYTRES,NREST,XB(1,1,IZ),ly1)
1442 100 CONTINUE
1443 CALL MXM (XB,NXY1,IZTRES,NZREST,Y,lz1)
1444 ENDIF
1445C
1446 RETURN
1447 END
1448C
1449 subroutine map31 (y,x,iel)
1450C---------------------------------------------------------------
1451C
1452C Map the elemental array X from mesh M3 to mesh M1
1453C
1454C---------------------------------------------------------------
1455C
1456 INCLUDE 'SIZE'
1457 INCLUDE 'GEOM'
1458 INCLUDE 'IXYZ'
1459 INCLUDE 'INPUT'
1460C
1461 REAL X(LX3,LY3,LZ3)
1462 REAL Y(LX1,LY1,LZ1)
1463C
1464 COMMON /CTMP0/ XA(LX1,LY3,LZ3), XB(LX1,LY1,LZ3)
1465C
1466 NYZ3 = ly3*lz3
1467 NXY1 = lx1*ly1
1468C
1469C Use the appropriate derivative- and interpolation operator in
1470C the y-direction (= radial direction if axisymmetric).
1471C
1472 IF (IFAXIS) THEN
1473 ly31 = ly1*ly3
1474 IF (IFRZER(IEL)) CALL COPY (IYTM31,IATM31,ly31)
1475 IF (.NOT.IFRZER(IEL)) CALL COPY (IYTM31,ICTM31,ly31)
1476 ENDIF
1477C
1478 IF (IF3D) THEN
1479 CALL MXM (IXM31,lx1,X,lx3,XA,NYZ3)
1480 DO 100 IZ=1,lz3
1481 CALL MXM (XA(1,1,IZ),lx1,IYTM31,ly3,XB(1,1,IZ),ly1)
1482 100 CONTINUE
1483 CALL MXM (XB,NXY1,IZTM31,lz3,Y,lz1)
1484 ELSE
1485 CALL MXM (IXM31,lx1,X,lx3,XA,NYZ3)
1486 CALL MXM (XA,lx1,IYTM31,ly3,Y,ly1)
1487 ENDIF
1488C
1489 RETURN
1490 END
1491C
1492 subroutine map13 (y,x,iel)
1493C---------------------------------------------------------------
1494C
1495C Map the elemental array X from mesh M1 to mesh M3
1496C
1497C---------------------------------------------------------------
1498C
1499 INCLUDE 'SIZE'
1500 INCLUDE 'GEOM'
1501 INCLUDE 'IXYZ'
1502 INCLUDE 'INPUT'
1503C
1504 REAL X(LX1,LY1,LZ1)
1505 REAL Y(LX3,LY3,LZ3)
1506C
1507 COMMON /CTMP0/ XA(LX3,LY1,LZ1), XB(LX3,LY3,LZ1)
1508C
1509 NYZ1 = ly1*lz1
1510 NXY3 = lx3*ly3
1511C
1512C Use the appropriate derivative- and interpolation operator in
1513C the y-direction (= radial direction if axisymmetric).
1514C
1515 IF (IFAXIS) THEN
1516 ly13 = ly1*ly3
1517 IF (IFRZER(IEL)) CALL COPY (IYTM13,IATM13,ly13)
1518 IF (.NOT.IFRZER(IEL)) CALL COPY (IYTM13,ICTM13,ly13)
1519 ENDIF
1520C
1521 CALL MXM (IXM13,lx3,X,lx1,XA,NYZ1)
1522 DO 100 IZ=1,lz1
1523 CALL MXM (XA(1,1,IZ),lx3,IYTM13,ly1,XB(1,1,IZ),ly3)
1524 100 CONTINUE
1525 CALL MXM (XB,NXY3,IZTM13,lz1,Y,lz3)
1526C
1527 RETURN
1528 END
1529 subroutine map12 (y,x,iel)
1530C---------------------------------------------------------------
1531C
1532C Map the elemental array X from mesh M1 to mesh M2
1533C
1534C---------------------------------------------------------------
1535C
1536 INCLUDE 'SIZE'
1537 INCLUDE 'GEOM'
1538 INCLUDE 'IXYZ'
1539 INCLUDE 'INPUT'
1540C
1541 REAL X(LX1,LY1,LZ1)
1542 REAL Y(LX2,LY2,LZ2)
1543C
1544 COMMON /CTMP00/ XA(LX2,LY1,LZ1), XB(LX2,LY2,LZ1)
1545C
1546 NYZ1 = ly1*lz1
1547 NXY2 = lx2*ly2
1548C
1549C Use the appropriate derivative- and interpolation operator in
1550C the y-direction (= radial direction if axisymmetric).
1551C
1552 IF (IFAXIS) THEN
1553 ly12 = ly1*ly2
1554 IF (IFRZER(IEL)) CALL COPY (IYTM12,IATM12,ly12)
1555 IF (.NOT.IFRZER(IEL)) CALL COPY (IYTM12,ICTM12,ly12)
1556 ENDIF
1557C
1558 CALL MXM (IXM12,lx2,X,lx1,XA,NYZ1)
1559 DO 100 IZ=1,lz1
1560 CALL MXM (XA(1,1,IZ),lx2,IYTM12,ly1,XB(1,1,IZ),ly2)
1561 100 CONTINUE
1562 CALL MXM (XB,NXY2,IZTM12,lz1,Y,lz2)
1563C
1564 RETURN
1565 END
1566C
1567 subroutine map21t (y,x,iel)
1568C---------------------------------------------------------------
1569C
1570C Map the elemental array X from mesh M2 to mesh M1 (Y)
1571C
1572C---------------------------------------------------------------
1573C
1574 INCLUDE 'SIZE'
1575 INCLUDE 'GEOM'
1576 INCLUDE 'IXYZ'
1577 INCLUDE 'INPUT'
1578C
1579 REAL X(LX2,LY2,LZ2)
1580 REAL Y(LX1,LY1,LZ1)
1581C
1582 COMMON /CTMP0/ XA(LX1,LY2,LZ2), XB(LX1,LY1,LZ2)
1583C
1584 NYZ2 = ly2*lz2
1585 NXY1 = lx1*ly1
1586 NXYZ = lx1*ly1*lz1
1587C
1588C Use the appropriate derivative- and interpolation operator in
1589C the y-direction (= radial direction if axisymmetric).
1590C
1591 IF (IFSPLIT) THEN
1592 CALL COPY(Y,X,NXYZ)
1593 RETURN
1594 ENDIF
1595C
1596 IF (IF3D) THEN
1597 CALL MXM (IXM21,lx1,X,lx2,XA,NYZ2)
1598 DO 100 IZ=1,lz2
1599 CALL MXM (XA(1,1,IZ),lx1,IYTM21,ly2,XB(1,1,IZ),ly1)
1600 100 CONTINUE
1601 CALL MXM (XB,NXY1,IZTM21,lz2,Y,lz1)
1602 ELSE
1603 CALL MXM (IXM21,lx1,X,lx2,XA,NYZ2)
1604 CALL MXM (XA,lx1,IYTM21,ly2,Y,ly1)
1605 ENDIF
1606 RETURN
1607 END
1608 subroutine map21e (y,x,iel)
1609C---------------------------------------------------------------
1610C
1611C Map the elemental array X from mesh M2 to mesh M1
1612C
1613C---------------------------------------------------------------
1614C
1615 INCLUDE 'SIZE'
1616 INCLUDE 'GEOM'
1617 INCLUDE 'IXYZ'
1618 INCLUDE 'INPUT'
1619C
1620 REAL X(LX2,LY2,LZ2)
1621 REAL Y(LX1,LY1,LZ1)
1622C
1623 COMMON /CTMP0/ XA(LX1,LY2,LZ2), XB(LX1,LY1,LZ2)
1624C
1625 NYZ2 = ly2*lz2
1626 NXY1 = lx1*ly1
1627C
1628C Use the appropriate derivative- and interpolation operator in
1629C the y-direction (= radial direction if axisymmetric).
1630C
1631 IF (IFAXIS) THEN
1632 ly21 = ly1*ly2
1633 IF (IFRZER(IEL)) CALL COPY (IYM12,IAM12,ly21)
1634 IF (.NOT.IFRZER(IEL)) CALL COPY (IYM12,ICM12,ly21)
1635 ENDIF
1636C
1637 CALL MXM (IXTM12,lx1,X,lx2,XA,NYZ2)
1638 DO 100 IZ=1,lz2
1639 CALL MXM (XA(1,1,IZ),lx1,IYM12,ly2,XB(1,1,IZ),ly1)
1640 100 CONTINUE
1641 CALL MXM (XB,NXY1,IZM12,lz2,Y,lz1)
1642C
1643 RETURN
1644 END
1645c-----------------------------------------------------------------------
1646 subroutine out_xyz_el(x,y,z,e)
1647 include 'SIZE'
1648 integer e
1649 real x(1),y(1),z(1)
1650c
1651 call out_fld_el(x,e,'XQ')
1652 call out_fld_el(y,e,'YQ')
1653 call out_fld_el(z,e,'ZQ')
1654c
1655 return
1656 end
1657c-----------------------------------------------------------------------
1658 subroutine out_fld_el(x,e,c2)
1659 include 'SIZE'
1660 real x(lx1,ly1,lz1,lelt)
1661 integer e
1662 character*2 c2
1663c
1664 write(6,1) c2,e
1665 nx8 = min(lx1,8)
1666 do k=1,lz1
1667 do j=1,ly1
1668 write(6,1) c2,e,(x(i,j,k,e),i=1,nx8)
1669 enddo
1670 enddo
1671 1 format(a2,i6,1p8e11.3)
1672 return
1673 end
1674c-----------------------------------------------------------------------
1675 subroutine outxm3j(xm3,ym3,jm3)
1676 include 'SIZE'
1677 include 'TOTAL'
1678
1679 real xm3(lx1,ly1,lz1,lelv)
1680 real ym3(lx1,ly1,lz1,lelv)
1681 real jm3(lx1,ly1,lz1,lelv)
1682
1683 integer e
1684
1685 do e=1,nelt
1686 write(6,*) e,nelfld(e),iftmsh(e),' iftmsh'
1687 call outmat(xm3(1,1,1,e),lx3,ly3,' xm3 ',e)
1688 call outmat(ym3(1,1,1,e),lx3,ly3,' ym3 ',e)
1689 call outmat(jm3(1,1,1,e),lx3,ly3,' jm3 ',e)
1690 enddo
1691
1692 return
1693 end
1694c-----------------------------------------------------------------------
1695 SUBROUTINE INVMT(A,B,AA,N)
1696C
1697 REAL A(N,N),AA(N,N),B(N,N)
1698 INTEGER INDX(100)
1699C
1700 NN = N*N
1701 DO 12 I=1,N
1702 DO 11 J=1,N
1703 B(I,J) = 0.0
1704 11 CONTINUE
1705 B(I,I) = 1.0
1706 12 CONTINUE
1707C
1708 CALL COPY (AA,A,NN)
1709 CALL LUDCMP(AA,N,N,INDX,D)
1710 DO 13 J=1,N
1711 CALL LUBKSB(AA,N,N,INDX,B(1,J))
1712 13 CONTINUE
1713C
1714 RETURN
1715 END
1716
1717 SUBROUTINE LUBKSB(A,N,NP,INDX,B)
1718 REAL A(NP,NP),B(N)
1719 INTEGER INDX(N)
1720 II=0
1721 DO 12 I=1,N
1722 LL=INDX(I)
1723 SUM=B(LL)
1724 B(LL)=B(I)
1725 IF (II.NE.0)THEN
1726 DO 11 J=II,I-1
1727 SUM=SUM-A(I,J)*B(J)
172811 CONTINUE
1729 ELSE IF (SUM.NE.0.0) THEN
1730 II=I
1731 ENDIF
1732 B(I)=SUM
173312 CONTINUE
1734 DO 14 I=N,1,-1
1735 SUM=B(I)
1736 IF(I.LT.N)THEN
1737 DO 13 J=I+1,N
1738 SUM=SUM-A(I,J)*B(J)
173913 CONTINUE
1740 ENDIF
1741 B(I)=SUM/A(I,I)
174214 CONTINUE
1743 RETURN
1744 END
1745 SUBROUTINE LUDCMP(A,N,NP,INDX,D)
1746 PARAMETER (NMAX=100,TINY=1.0E-20)
1747 REAL A(NP,NP),VV(NMAX)
1748 INTEGER INDX(N)
1749 D=1.0
1750 DO 12 I=1,N
1751 AAMAX=0.0
1752 DO 11 J=1,N
1753 IF (ABS(A(I,J)).GT.AAMAX) AAMAX=ABS(A(I,J))
175411 CONTINUE
1755 IF (AAMAX.EQ.0.0) THEN
1756 write(6,*) 'Singular matrix.'
1757 call exitt
1758 ENDIF
1759 VV(I)=1.0/AAMAX
176012 CONTINUE
1761 DO 19 J=1,N
1762 IF (J.GT.1) THEN
1763 DO 14 I=1,J-1
1764 SUM=A(I,J)
1765 IF (I.GT.1)THEN
1766 DO 13 K=1,I-1
1767 SUM=SUM-A(I,K)*A(K,J)
176813 CONTINUE
1769 A(I,J)=SUM
1770 ENDIF
177114 CONTINUE
1772 ENDIF
1773 AAMAX=0.0
1774 DO 16 I=J,N
1775 SUM=A(I,J)
1776 IF (J.GT.1)THEN
1777 DO 15 K=1,J-1
1778 SUM=SUM-A(I,K)*A(K,J)
177915 CONTINUE
1780 A(I,J)=SUM
1781 ENDIF
1782 DUM=VV(I)*ABS(SUM)
1783 IF (DUM.GE.AAMAX) THEN
1784 IMAX=I
1785 AAMAX=DUM
1786 ENDIF
178716 CONTINUE
1788 IF (J.NE.IMAX)THEN
1789 DO 17 K=1,N
1790 DUM=A(IMAX,K)
1791 A(IMAX,K)=A(J,K)
1792 A(J,K)=DUM
179317 CONTINUE
1794 D=-D
1795 VV(IMAX)=VV(J)
1796 ENDIF
1797 INDX(J)=IMAX
1798 IF(J.NE.N)THEN
1799 IF(A(J,J).EQ.0.)A(J,J)=TINY
1800 DUM=1.0/A(J,J)
1801 DO 18 I=J+1,N
1802 A(I,J)=A(I,J)*DUM
180318 CONTINUE
1804 ENDIF
180519 CONTINUE
1806 IF(A(N,N).EQ.0.0)A(N,N)=TINY
1807 RETURN
1808 END
1809c-----------------------------------------------------------------------
1810 subroutine set_unr
1811 include 'SIZE'
1812 include 'INPUT'
1813 include 'GEOM'
1814 include 'TOPOL'
1815
1816 integer e,f,pf
1817
1818 nface = 2*ldim
1819 call dsset(lx1,ly1,lz1)
1820
1821 do e=1,nelt
1822 do f=1,nface
1823
1824 pf = eface1(f)
1825 js1 = skpdat(1,pf)
1826 jf1 = skpdat(2,pf)
1827 jskip1 = skpdat(3,pf)
1828 js2 = skpdat(4,pf)
1829 jf2 = skpdat(5,pf)
1830 jskip2 = skpdat(6,pf)
1831
1832 i = 0
1833 if (if3d) then
1834 do j2=js2,jf2,jskip2
1835 do j1=js1,jf1,jskip1
1836 i = i+1
1837 a = area(i,1,f,e)/jacm1(j1,j2,1,e)
1838 unr(i,f,e) = a * ( rxm1(j1,j2,1,e)*unx(i,1,f,e)
1839 $ + rym1(j1,j2,1,e)*uny(i,1,f,e)
1840 $ + rzm1(j1,j2,1,e)*unz(i,1,f,e) )
1841 uns(i,f,e) = a * ( sxm1(j1,j2,1,e)*unx(i,1,f,e)
1842 $ + sym1(j1,j2,1,e)*uny(i,1,f,e)
1843 $ + szm1(j1,j2,1,e)*unz(i,1,f,e) )
1844 unt(i,f,e) = a * ( txm1(j1,j2,1,e)*unx(i,1,f,e)
1845 $ + tym1(j1,j2,1,e)*uny(i,1,f,e)
1846 $ + tzm1(j1,j2,1,e)*unz(i,1,f,e) )
1847 enddo
1848 enddo
1849 else
1850 do j2=js2,jf2,jskip2
1851 do j1=js1,jf1,jskip1
1852 i = i+1
1853 a = area(i,1,f,e)/jacm1(j1,j2,1,e)
1854 unr(i,f,e) = a * ( rxm1(j1,j2,1,e)*unx(i,1,f,e)
1855 $ + rym1(j1,j2,1,e)*uny(i,1,f,e) )
1856 uns(i,f,e) = a * ( sxm1(j1,j2,1,e)*unx(i,1,f,e)
1857 $ + sym1(j1,j2,1,e)*uny(i,1,f,e) )
1858 unt(i,f,e) = 0.
1859 enddo
1860 enddo
1861 endif
1862 enddo
1863 enddo
1864
1865 return
1866 end
1867c-----------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.