OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
upd_glob_k.F
Go to the documentation of this file.
1Copyright> OpenRadioss
2Copyright> Copyright (C) 1986-2025 Altair Engineering Inc.
3Copyright>
4Copyright> This program is free software: you can redistribute it and/or modify
5Copyright> it under the terms of the GNU Affero General Public License as published by
6Copyright> the Free Software Foundation, either version 3 of the License, or
7Copyright> (at your option) any later version.
8Copyright>
9Copyright> This program is distributed in the hope that it will be useful,
10Copyright> but WITHOUT ANY WARRANTY; without even the implied warranty of
11Copyright> MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12Copyright> GNU Affero General Public License for more details.
13Copyright>
14Copyright> You should have received a copy of the GNU Affero General Public License
15Copyright> along with this program. If not, see <https://www.gnu.org/licenses/>.
16Copyright>
17Copyright>
18Copyright> Commercial Alternative: Altair Radioss Software
19Copyright>
20Copyright> As an alternative to this open-source version, Altair also offers Altair Radioss
21Copyright> software under a commercial license. Contact Altair to discuss further if the
22Copyright> commercial version may interest you: https://www.altair.com/radioss/.
23!||====================================================================
24!|| upd_glob_k ../engine/source/implicit/upd_glob_k.F
25!||--- called by ------------------------------------------------------
26!|| imp_buck ../engine/source/implicit/imp_buck.F
27!|| imp_chkm ../engine/source/implicit/imp_solv.f
28!|| imp_k_eig ../engine/stub/imp_k_eig.F
29!|| imp_solv ../engine/source/implicit/imp_solv.F
30!||--- calls -----------------------------------------------------
31!|| bc_imp0 ../engine/source/constraints/general/bcs/bc_imp0.F
32!|| bc_imp1 ../engine/source/constraints/general/bcs/bc_imp0.F
33!|| bc_impa ../engine/source/constraints/general/bcs/bc_imp0.F
34!|| condens_b ../engine/source/implicit/upd_glob_k.F
35!|| condens_k ../engine/source/implicit/upd_glob_k.F
36!|| dim_fvbcl ../engine/source/constraints/general/impvel/fv_imp0.F
37!|| fv_imp0 ../engine/source/constraints/general/impvel/fv_imp0.F
38!|| fv_impl ../engine/source/constraints/general/impvel/fv_imp0.F
39!|| fv_rw0 ../engine/source/constraints/general/impvel/fv_imp0.F
40!|| fvbc_allo ../engine/source/constraints/general/impvel/fv_imp0.F
41!|| fvbc_impl ../engine/source/constraints/general/impvel/fv_imp0.F
42!|| i2_imp0 ../engine/source/interfaces/interf/i2_imp0.F
43!|| ini_dofspc ../engine/source/implicit/upd_glob_k.F
44!|| rbe2_imp0 ../engine/source/constraints/general/rbe2/rbe2_imp0.f
45!|| rbe3_imp0 ../engine/source/constraints/general/rbe3/rbe3_imp0.F
46!|| rby_imp0 ../engine/source/constraints/general/rbody/rby_imp0.F
47!|| rm_imp0 ../engine/source/model/remesh/rm_imp0.F
48!|| upd_aspc ../engine/source/constraints/general/bcs/bc_imp0.F
49!||--- uses -----------------------------------------------------
50!|| imp_fvbcl ../engine/share/modules/impbufdef_mod.f
51!|| intbufdef_mod ../common_source/modules/interfaces/intbufdef_mod.F90
52!||====================================================================
53 SUBROUTINE upd_glob_k(ICODT ,ICODR ,ISKEW ,IBFV ,NPC ,
54 1 TF ,VEL ,XFRAME ,
55 2 RBY ,X ,SKEW ,LPBY ,NPBY ,
56 3 ITAB ,WEIGHT,MS ,IN ,NRBYAC ,
57 4 IRBYAC,NSC ,ISIJ ,NMC ,IMIJ ,
58 5 NSS ,ISS ,NINT2 ,IINT2 ,NSC2 ,
59 6 ISIJ2 ,NSS2 ,ISS2 ,IPARI ,INTBUF_TAB,
60 7 NDDL ,NNZ ,IADK ,JDIK ,
61 8 DIAG_K,LT_K ,NDOF ,IDDL ,IKC ,
62 9 UD ,B ,NKUD ,IKUD ,BKUD ,
63 A NMC2 ,IMIJ2 ,NT_RW ,RD ,LJ ,
64 B IRBE3 ,LRBE3 ,FRBE3 ,ISS3 ,IRBE2 ,
65 C LRBE2 ,ISB2 ,NSRB2 )
66C-----------------------------------------------
67C M o d u l e s
68C-----------------------------------------------
69 USE imp_fvbcl
70 USE intbufdef_mod
71C-----------------------------------------------
72C I m p l i c i t T y p e s
73C-----------------------------------------------
74#include "implicit_f.inc"
75C-----------------------------------------------
76C C o m m o n B l o c k s
77C-----------------------------------------------
78#include "com04_c.inc"
79#include "param_c.inc"
80#include "remesh_c.inc"
81#include "impl1_c.inc"
82C-----------------------------------------------
83C D u m m y A r g u m e n t s
84C-----------------------------------------------
85 INTEGER NPC(*),IBFV(NIFV,*),
86 . ICODT(*),ICODR(*),ISKEW(*),NSC(*),ISIJ(*),
87 . NMC,IMIJ(*),NSS(*),ISS(*),NINT2 ,IINT2(*),NT_RW,
88 . NSC2(*),ISIJ2(*),NSS2(*),ISS2(*),NKUD(*) ,IKUD(*),
89 . NMC2,IMIJ2(*),LJ(*)
90 INTEGER WEIGHT(*),LPBY(*),NPBY(NNPBY,*),ITAB(*),
91 . IPARI(NPARI,*), NRBYAC,IRBYAC(*),
92 . NDDL,NNZ,IADK(*),JDIK(*),NDOF(*),IDDL(*),IKC(*),
93 . IRBE3(*) ,LRBE3(*),ISS3(*),IRBE2(*),LRBE2(*),
94 . ISB2(*),NSRB2(*)
95 my_real
96 . RBY(NRBY,*) ,X(3,*) ,SKEW(*),IN(*),MS(*)
97 my_real
98 . tf(*),vel(lfxvelr,*),ud(3,*),diag_k(*),lt_k(*),
99 . b(*) ,bkud(*),rd(3,*),xframe(nxframe,*),frbe3(*)
100
101 TYPE(intbuf_struct_) INTBUF_TAB(*)
102C-----------------------------------------------
103C L o c a l V a r i a b l e s
104C-----------------------------------------------
105 INTEGER I,J,K
106 my_real,
107 . DIMENSION(:),ALLOCATABLE :: dofspc
108C----IKC(NDDL) : 1->fix_global,2->Ud_global,3,4->fv_rw,5->Int2_0,6->Int2_1,7->RB-
109C----a faire (8->fix_local,9->Ud_local)-----10,11->sliding rw---12->remesh,13->RBE3,14->AUTO-SPC
110C----15->AUTO-SPC-in local sys,16->RBE2,
111C----No ddl -> IDDL(NDDL)+NDOF(NUMNOD)-MIN(IKC,1)---
112 IF(nint2>0)THEN
113 CALL i2_imp0(nint2 ,iint2 ,
114 1 ipari,intbuf_tab,x,ms ,in ,
115 1 nmc2 ,imij2 ,itab ,
116 2 nsc2 ,isij2 ,nss2 ,iss2,
117 3 weight,ikc ,ndof ,nddl,iddl ,
118 4 iadk ,jdik ,diag_k ,lt_k ,b )
119 ENDIF
120c
121 IF(nrbe2>0)THEN
122 CALL rbe2_imp0(
123 1 irbe2 ,lrbe2 ,x ,nsrb2 ,isb2 ,
124 2 ikc ,ndof ,iddl ,iadk ,jdik ,
125 3 diag_k ,lt_k ,b ,weight ,itab ,
126 4 skew )
127 ENDIF
128 IF(nrbe3>0)THEN
129 CALL rbe3_imp0(
130 1 irbe3 ,lrbe3 ,frbe3 ,x ,skew ,
131 2 iss3 ,ikc ,ndof ,iddl ,iadk ,
132 3 jdik ,diag_k,lt_k ,b ,weight ,
133 4 itab )
134 ENDIF
135 IF(nrbyac>0)THEN
136 CALL rby_imp0(x ,rby ,lpby ,npby ,skew ,
137 1 nrbyac,irbyac,nsc ,isij ,nmc ,
138 2 imij ,nss ,iss ,iskew ,itab ,
139 3 weight,ms ,in ,
140 4 nddl ,iadk ,jdik ,diag_k,
141 5 lt_k ,ndof ,iddl ,ikc ,b )
142 ENDIF
143 IF(nadmesh/=0)THEN
144 CALL rm_imp0(
145 4 nddl ,iadk ,jdik ,diag_k ,lt_k ,
146 5 ndof ,iddl ,ikc ,b ,itab )
147 ENDIF
148C -----enforce b.c. FV va tenir compte b.c.---------------------
149 CALL bc_imp0(icodt ,icodr,iskew,ikc,ndof,iddl)
150C-----case coupling FV+BCS using SKEW ----
151 nfvbcl=0
152 IF(nfxvel > 0) THEN
153 CALL dim_fvbcl(ibfv ,lj ,iskew ,icodt ,icodr ,
154 1 nddl ,iddl ,ikc ,iadk ,jdik ,
155 2 skew ,nfvbcl,nkud_l )
156c NFVBCL=0
157 IF(nfvbcl > 0 ) THEN
158 CALL fvbc_allo()
159 CALL fvbc_impl(ibfv ,skew ,xframe ,lj ,iddl ,
160 1 ikc ,ndof ,iadk ,jdik ,diag_k ,
161 2 lt_k ,ud ,rd ,b ,nddl ,
162 3 icodt ,icodr ,ict_1 ,icr_1,nkud_1,
163 4 ikud_1 ,bkud_1 )
164 END if!(NFVBCL > 0) THEN
165 ENDIF
166C
167 IF(nfvbcl > 0) THEN
168 CALL bc_imp1(ict_1 ,icr_1 ,iskew ,skew ,ikc ,
169 1 ndof ,iddl ,iadk ,jdik ,diag_k,
170 2 lt_k )
171 ELSE
172 CALL bc_imp1(icodt ,icodr ,iskew ,skew ,ikc ,
173 1 ndof ,iddl ,iadk ,jdik ,diag_k,
174 2 lt_k )
175 END IF
176 IF(nfxvel>0) THEN
177 CALL fv_impl(ibfv ,skew ,xframe ,lj ,iddl ,
178 1 ikc ,ndof ,iadk ,jdik ,diag_k ,
179 2 lt_k ,ud ,rd ,b )
180 CALL fv_imp0(iddl ,ikc ,ndof ,iadk ,jdik ,
181 1 diag_k ,lt_k ,ud ,nkud ,ikud ,
182 2 bkud ,nddl ,rd )
183 ENDIF
184C
185 IF(nt_rw>0) THEN
186 CALL fv_rw0(iddl ,ikc ,ndof ,iadk ,jdik ,
187 1 diag_k ,lt_k ,ud ,b )
188 ENDIF
189 IF(iautspc>0) THEN
190 ALLOCATE(dofspc(numnod))
191 CALL ini_dofspc(
192 1 npby ,lpby ,nrbyac ,irbyac ,nint2 ,
193 2 iint2 ,ipari ,intbuf_tab,ndof ,irbe3 ,
194 3 lrbe3 ,irbe2 ,lrbe2 ,x ,dofspc )
195 DO i=1,nddl
196 IF (ikc(i)==14.OR.ikc(i)==15) ikc(i)=0
197 END DO
198 CALL upd_aspc(nddl ,dofspc,iddl ,ikc ,itab ,
199 . iadk ,jdik ,diag_k,lt_k )
200 CALL bc_impa(iadk ,jdik ,diag_k,lt_k ,ndof ,
201 1 iddl ,ikc )
202 DEALLOCATE(dofspc)
203 ENDIF
204C -----update---condense-IKC>=1------------------------------------
205 CALL condens_b(nddl ,ikc ,diag_k)
206 CALL condens_k(nddl ,nnz ,iadk ,jdik ,lt_k ,
207 1 ikc )
208C
209 RETURN
210 END
211!||====================================================================
212!|| condens_k ../engine/source/implicit/upd_glob_k.F
213!||--- called by ------------------------------------------------------
214!|| upd_glob_k ../engine/source/implicit/upd_glob_k.F
215!|| upd_int_k ../engine/source/implicit/upd_glob_k.F
216!||--- calls -----------------------------------------------------
217!|| condens_k0 ../engine/source/implicit/upd_glob_k.F
218!||====================================================================
219 SUBROUTINE condens_k(
220 1 NDDL ,NNZ ,IADK ,JDIK ,
221 2 LT_K ,IKC )
222C-----------------------------------------------
223C I m p l i c i t T y p e s
224C-----------------------------------------------
225#include "implicit_f.inc"
226C-----------------------------------------------
227C C o m m o n B l o c k s
228C-----------------------------------------------
229#include "impl1_c.inc"
230C-----------------------------------------------
231C D u m m y A r g u m e n t s
232C-----------------------------------------------
233 INTEGER NDDL,NNZ,IADK(*),JDIK(*),IKC(*)
234 my_real
235 . LT_K(*)
236C-----------------------------------------------
237C L o c a l V a r i a b l e s
238C-----------------------------------------------
239 INTEGER I,J,K,NL,ND,IFIX(NDDL),NFT,JLT,NZF,NZL
240C -----update---condense-IKC>=1------------------------------------
241 ND = nddl
242 nddl=0
243 nnz=0
244 nl = 0
245C
246 DO nft = 0 , nd-1 , nnsiz
247 jlt = min( nnsiz, nd - nft )
248 nzf=iadk(nft+1)-1
249 nzl=iadk(jlt+nft+1)-iadk(nft+1)
250 CALL condens_k0(
251 1 nddl ,nnz ,iadk ,jdik ,nl ,
252 2 lt_k ,ikc ,ifix ,nft ,jlt ,
253 3 nzf ,nzl )
254 ENDDO
255 iadk(nddl+1)=nnz+1
256C
257 DO j = 1,nnz
258 k= jdik(j)
259 jdik(j) = jdik(j) -ifix(k)
260 ENDDO
261C
262 RETURN
263 END
264!||====================================================================
265!|| condens_k0 ../engine/source/implicit/upd_glob_k.F
266!||--- called by ------------------------------------------------------
267!|| condens_k ../engine/source/implicit/upd_glob_k.F
268!||====================================================================
269 SUBROUTINE condens_k0(
270 1 NDDL ,NNZ ,IADK ,JDIK ,NL ,
271 2 LT_K ,IKC ,IFIX ,NDF ,NDL ,
272 3 NZF ,NZL )
273C-----------------------------------------------
274C I m p l i c i t T y p e s
275C-----------------------------------------------
276#include "implicit_f.inc"
277C-----------------------------------------------
278C D u m m y A r g u m e n t s
279C-----------------------------------------------
280 INTEGER NDDL,NNZ,IADK(*),JDIK(*),IKC(*),IFIX(*),NL,
281 . NDF,NDL,NZF,NZL
282 my_real
283 . LT_K(*)
284C-----------------------------------------------
285C L o c a l V a r i a b l e s
286C-----------------------------------------------
287 INTEGER IADK1(NDL+1),JDIK1(NZL)
288 INTEGER I,J,K,I1,J1,ND,NI,JI
289 my_real
290 . LT_K1(NZL)
291C -----update---condense-IKC>=1------------------------------------
292 DO NI = 1,ndl+1
293 i=ni+ndf
294 iadk1(ni)=iadk(i)
295 ENDDO
296 DO ni = 1,nzl
297 i=ni+nzf
298 jdik1(ni)=jdik(i)
299 lt_k1(ni)=lt_k(i)
300 ENDDO
301C
302 DO ni = 1,ndl
303 i=ni+ndf
304 ifix(i) = nl
305 IF (ikc(i)<1) THEN
306 nddl=nddl+1
307 iadk(nddl)=nnz+1
308 DO j = iadk1(ni),iadk1(ni+1)-1
309 ji=j-nzf
310 k= jdik1(ji)
311 IF (ikc(k)<1) THEN
312 nnz=nnz+1
313 jdik(nnz)=k
314 lt_k(nnz)=lt_k1(ji)
315 ENDIF
316 ENDDO
317 ELSE
318 nl=nl+1
319 ENDIF
320 ENDDO
321C
322 RETURN
323 END
324!||====================================================================
325!|| condens_ind ../engine/source/implicit/upd_glob_k.F
326!||--- called by ------------------------------------------------------
327!|| upd_fr_k ../engine/source/mpi/implicit/imp_fri.F
328!||====================================================================
329 SUBROUTINE condens_ind(NDDL ,NNZ ,IADK ,JDIK ,IKC )
330C-----------------------------------------------
331C I m p l i c i t T y p e s
332C-----------------------------------------------
333#include "implicit_f.inc"
334C-----------------------------------------------
335C D u m m y A r g u m e n t s
336C-----------------------------------------------
337 INTEGER NDDL,NNZ,IADK(*),JDIK(*),IKC(*)
338C-----------------------------------------------
339C L o c a l V a r i a b l e s
340C-----------------------------------------------
341 INTEGER IADK1(NDDL+1),JDIK1(NNZ),IFIX(NDDL)
342 INTEGER I,J,K,I1,J1,ND,NI,JI,NL
343C -----update---condense-IKC>=1------------------------------------
344 IF (NDDL==0) return
345 DO i = 1,nddl+1
346 iadk1(i)=iadk(i)
347 ENDDO
348 DO ni = 1,nnz
349 jdik1(ni)=jdik(ni)
350 ENDDO
351C
352 nd = 0
353 nnz = 0
354 nl = 0
355 DO i = 1,nddl
356 ifix(i) = nl
357 IF (ikc(i)<1) THEN
358 nd=nd+1
359 iadk(nd)=nnz+1
360 DO j = iadk1(i),iadk1(i+1)-1
361 k= jdik1(j)
362 IF (ikc(k)<1) THEN
363 nnz=nnz+1
364 jdik(nnz)=k
365 ENDIF
366 ENDDO
367 ELSE
368 nl=nl+1
369 ENDIF
370 ENDDO
371 iadk(nd+1)=nnz+1
372 nddl=nd
373C
374 DO j = 1,nnz
375 k= jdik(j)
376 jdik(j) = jdik(j) -ifix(k)
377 ENDDO
378C
379 RETURN
380 END
381!||====================================================================
382!|| condens_b ../engine/source/implicit/upd_glob_k.F
383!||--- called by ------------------------------------------------------
384!|| d_to_u ../engine/source/implicit/produt_v.F
385!|| ext_rhs ../engine/source/implicit/upd_glob_k.F
386!|| get_fext ../engine/source/implicit/imp_solv.F
387!|| imp_fhht1 ../engine/source/implicit/imp_dyna.F
388!|| imp_solv ../engine/source/implicit/imp_solv.F
389!|| produt_u ../engine/source/implicit/produt_v.F
390!|| produt_u2 ../engine/source/implicit/produt_v.F
391!|| produt_uh ../engine/source/implicit/produt_v.F
392!|| produt_uh2 ../engine/source/implicit/produt_v.F
393!|| produt_uhp ../engine/source/implicit/produt_v.F
394!|| produt_uhp2 ../engine/source/implicit/produt_v.f
395!|| upd_glob_k ../engine/source/implicit/upd_glob_k.F
396!|| upd_int_k ../engine/source/implicit/upd_glob_k.f
397!|| upd_rhs ../engine/source/implicit/upd_glob_k.F
398!||====================================================================
399 SUBROUTINE condens_b(NDDL ,IKC ,B )
400C-----------------------------------------------
401C I m p l i c i t T y p e s
402C-----------------------------------------------
403#include "implicit_f.inc"
404C-----------------------------------------------
405C D u m m y A r g u m e n t s
406C-----------------------------------------------
407 INTEGER NDDL,IKC(*)
408 my_real
409 . B(*)
410C-----------------------------------------------
411C L o c a l V a r i a b l e s
412C-----------------------------------------------
413 INTEGER I,ND
414 my_real
415 . B1(NDDL)
416C -----update---condense-IKC>=1------------------------------------
417 DO i = 1,nddl
418 b1(i)=b(i)
419 ENDDO
420C
421 nd=0
422C
423 DO i = 1,nddl
424 IF (ikc(i)<1) THEN
425 nd=nd+1
426 b(nd)=b1(i)
427 ENDIF
428 ENDDO
429C
430 RETURN
431 END
432!||====================================================================
433!|| upd_int_k ../engine/source/implicit/upd_glob_k.F
434!||--- called by ------------------------------------------------------
435!|| imp_int_k ../engine/source/implicit/imp_int_k.F
436!||--- calls -----------------------------------------------------
437!|| bc_imp0 ../engine/source/constraints/general/bcs/bc_imp0.F
438!|| bc_imp1 ../engine/source/constraints/general/bcs/bc_imp0.f
439!|| bc_impa ../engine/source/constraints/general/bcs/bc_imp0.F
440!|| condens_b ../engine/source/implicit/upd_glob_k.F
441!|| condens_k ../engine/source/implicit/upd_glob_k.F
442!|| cp_real ../engine/source/implicit/produt_v.f
443!|| fv_impi ../engine/source/constraints/general/impvel/fv_imp0.F
444!|| fv_impl ../engine/source/constraints/general/impvel/fv_imp0.F
445!|| fv_rw0 ../engine/source/constraints/general/impvel/fv_imp0.F
446!|| i2_impi ../engine/source/interfaces/interf/i2_imp0.F
447!|| rbe2_impi ../engine/source/constraints/general/rbe2/rbe2_imp0.F
448!|| rbe3_impi ../engine/source/constraints/general/rbe3/rbe3_imp0.F
449!|| rby_impi ../engine/source/constraints/general/rbody/rby_imp0.F
450!||--- uses -----------------------------------------------------
451!|| intbufdef_mod ../common_source/modules/interfaces/intbufdef_mod.F90
452!||====================================================================
453 SUBROUTINE upd_int_k(ICODT ,ICODR ,ISKEW ,IBFV ,NPC ,
454 1 TF ,VEL ,XFRAME ,
455 2 RBY ,X ,SKEW ,LPBY ,NPBY ,
456 3 ITAB ,WEIGHT,MS ,IN ,NRBYAC,
457 4 IRBYAC,NSS ,ISS ,IPARI ,INTBUF_TAB,
458 5 NINT2 ,IINT2 ,IAINT2 ,NSS2 ,
459 5 ISS2 ,NDDLI ,NNZI ,IADI ,JDII ,
460 6 DIAG_I ,LT_I ,IDDLI ,NDDL ,IADK ,
461 7 JDIK ,IKC ,DIAG_K,LT_K ,IDDL ,
462 8 NDOFI ,ITOK ,UD ,LB ,LUJ ,
463 9 NT_RW ,IRBE3 ,LRBE3 ,FRBE3 ,NSS3 ,
464 A ISS3 ,IRBE2 ,LRBE2 ,NSB2 ,ISB2 )
465C-----------------------------------------------
466C M o d u l e s
467C-----------------------------------------------
468 USE intbufdef_mod
469C-----------------------------------------------
470C I m p l i c i t T y p e s
471C-----------------------------------------------
472#include "implicit_f.inc"
473C-----------------------------------------------
474C C o m m o n B l o c k s
475C-----------------------------------------------
476#include "com04_c.inc"
477#include "param_c.inc"
478#include "impl1_c.inc"
479C-----------------------------------------------
480C D u m m y A r g u m e n t s
481C-----------------------------------------------
482 INTEGER NPC(*),IBFV(NIFV,*),LUJ(*),
483 . ICODT(*),ICODR(*),ISKEW(*),ITOK(*),NDDL,NT_RW
484 INTEGER WEIGHT(*),LPBY(*),NPBY(NNPBY,*),ITAB(*),
485 . IPARI(NPARI,*), NRBYAC,IRBYAC(*),
486 . IDDL(*),IKC(*),NSS(*),ISS(*),NSS2(*),ISS2(*),
487 . IADK(*),JDIK(*),NDDLI,NNZI,IADI(*),JDII(*),
488 . IDDLI(*),NDOFI(*),NINT2 ,IINT2(*),IAINT2(*)
489 INTEGER IRBE3(NRBE3L,*),LRBE3(*),NSS3(*),ISS3(*),
490 . IRBE2(*),LRBE2(*),NSB2(*),ISB2(*)
491 my_real
492 . RBY(NRBY,*) ,X(3,*) ,SKEW(*),IN(*),MS(*)
493 my_real
494 . TF(*), VEL(LFXVELR,*),DIAG_K(*),LT_K(*),
495 . DIAG_I(*),LT_I(*),LB(*),UD(3,*),XFRAME(NXFRAME,*),FRBE3(*)
496
497 TYPE(intbuf_struct_) INTBUF_TAB(*)
498C-----------------------------------------------
499C L o c a l V a r i a b l e s
500C-----------------------------------------------
501 INTEGER LUJI(NFXVEL),IBID
502 INTEGER I,J,N, IAD,NTY,NDDL0,IIC(NDDLI),ND,
503 . ITOK1(NDDLI),IDEN(NDDL),NNZK0,L1,L2,K,IKP,NN
504C REAL
505 my_real
506 . LBI(NDDLI),RBID,LB0(NDDLI)
507C--------diak=0->IIC=10--------------------------------------
508C
509 DO I=1,nddli
510 iad=itok(i)
511 iic(i)=0
512 IF (ikc(iad)==2.OR.ikc(iad)==9) iic(i)=ikc(iad)
513 IF (ikc(iad)==3.OR.ikc(iad)==4.OR.
514 . ikc(iad)==10.OR.ikc(iad)==11) iic(i)=ikc(iad)
515 IF (ikc(iad)==14.OR.ikc(iad)==15) iic(i)=ikc(iad)
516 IF (ikc(iad)<1) THEN
517 lbi(i)=lb(iad)
518 ELSE
519 lbi(i)=zero
520 ENDIF
521 itok1(i)=iad
522 ENDDO
523 IF(nint2>0.OR.nrbyac>0.OR.nrbe3>0)
524 . CALL cp_real(nddli,lbi,lb0)
525 nddl0=nddli
526 nnzk0=nnzi
527C -----update---due to kinematic conditions-----------------------------------
528C----IKC(NDDL) : 1->fix_global,2->Ud_global,3->fix_local,4->Ud_local,5->Int2_0,6->Int2_1,7->RB------
529C-----l'interface 2 ------
530 IF(nint2>0)THEN
531 CALL i2_impi(nint2 ,iint2 ,ipari ,intbuf_tab ,
532 1 x ,ms ,in ,nss2 ,iss2 ,
533 2 weight ,iic ,ndofi ,nddli ,iddli ,
534 3 iadi ,jdii ,diag_i,lt_i ,iaint2,
535 4 lb0 ,itab )
536 ENDIF
537 IF(nrbe2>0)THEN
538 CALL rbe2_impi(
539 1 irbe2 ,lrbe2 ,x ,skew ,
540 2 nsb2 ,isb2 ,iic ,ndofi ,iddli ,
541 3 iadi ,jdii ,diag_i ,lt_i ,lb0 ,
542 4 weight ,itab )
543 ENDIF
544C-----rbe3 ------
545 IF(nrbe3>0)THEN
546 CALL rbe3_impi(
547 1 irbe3 ,lrbe3 ,frbe3 ,x ,skew ,
548 2 nss3 ,iss3 ,iic ,ndofi ,iddli ,
549 3 iadi ,jdii ,diag_i,lt_i ,lb0 ,
550 4 weight ,itab )
551 ENDIF
552C-----RIGID BODIES-------
553 IF(nrbyac>0)THEN
554 CALL rby_impi(x ,rby ,lpby ,npby ,skew ,
555 1 nrbyac,irbyac,nss ,iss ,iskew ,
556 2 itab ,weight,ms ,in ,
557 3 nddli ,iadi ,jdii ,diag_i,
558 4 lt_i ,ndofi ,iddli ,iic ,lb0 )
559 ENDIF
560C------------
561 CALL bc_imp0(icodt ,icodr,iskew,iic,ndofi,iddli)
562 CALL bc_imp1(icodt ,icodr ,iskew ,skew ,iic ,
563 1 ndofi ,iddli ,iadi ,jdii ,diag_i,
564 2 lt_i )
565 IF(iautspc>0) THEN
566 CALL bc_impa(iadi ,jdii ,diag_i,lt_i ,ndofi ,
567 1 iddli ,iic )
568 ENDIF
569 IF(nfxvel>0) THEN
570 DO i=1,nfxvel
571 n=iabs(ibfv(1,i))
572 IF (ndofi(n)>0.AND.luj(i)<=3) THEN
573 luji(i) = luj(i)
574 ELSE
575 luji(i) = 0
576 ENDIF
577 ENDDO
578 CALL fv_impl(ibfv ,skew ,xframe ,luji ,iddli ,
579 1 iic ,ndofi ,iadi ,jdii ,diag_i ,
580 2 lt_i ,ud ,rbid ,lbi )
581 CALL fv_impi(iddli ,iic ,ndofi ,iadi ,jdii ,
582 1 diag_i ,lt_i ,ud ,lbi ,nddli )
583 ENDIF
584 IF(nt_rw>0) THEN
585 CALL fv_rw0(iddli ,iic ,ndofi ,iadi ,jdii ,
586 1 diag_i ,lt_i ,ud ,lbi )
587 ENDIF
588C -----diak=0------------------------------------
589 IF (intp_c>=0) THEN
590 DO i=1,nddli
591 IF (iic(i)<1.AND.diag_i(i)<=em20) iic(i)=10
592 ENDDO
593 ENDIF
594C -----update---LB------------------------------------
595 DO i=1,nddli
596 iad=itok(i)
597 IF (iic(i)<1) lb(iad)=lbi(i)
598 ENDDO
599C -----update---condense-IKC>=1------------------------------------
600 CALL condens_b(nddli ,iic ,diag_i)
601 CALL condens_k(nddli ,nnzi ,iadi ,jdii ,lt_i ,iic )
602C -----update---ITOK------------------------------------
603 nd=0
604 DO i = 1,nddl
605 IF (ikc(i)<1) THEN
606 nd=nd+1
607 iden(i)=nd
608 ELSE
609 iden(i)=-ikc(i)
610 ENDIF
611 ENDDO
612C
613 nd=0
614 DO i = 1,nddl0
615 IF (iic(i)<1) THEN
616 nd=nd+1
617 l1=itok1(i)
618 l2=iden(l1)
619 IF (l2<=0)WRITE(*,*)'---ERROR ITOK:---',i,l1,l2
620 itok(nd)=l2
621 ENDIF
622 ENDDO
623C
624 RETURN
625 END
626!||====================================================================
627!|| upd_rhs ../engine/source/implicit/upd_glob_k.F
628!||--- called by ------------------------------------------------------
629!|| imp_chkm ../engine/source/implicit/imp_solv.F
630!|| imp_solv ../engine/source/implicit/imp_solv.F
631!||--- calls -----------------------------------------------------
632!|| bc_impr1 ../engine/source/constraints/general/bcs/bc_imp0.F
633!|| condens_b ../engine/source/implicit/upd_glob_k.F
634!|| ext_rhs ../engine/source/implicit/upd_glob_k.F
635!|| fv_imprl ../engine/source/constraints/general/impvel/fv_imp0.F
636!|| fv_rwlr0 ../engine/source/constraints/general/rwall/srw_imp.F
637!|| i2_impr1 ../engine/source/interfaces/interf/i2_imp1.F
638!|| i2_impr2 ../engine/source/interfaces/interf/i2_imp1.F
639!|| imp_dycrb ../engine/source/implicit/imp_dyna.F
640!|| imp_dynar ../engine/source/implicit/imp_dyna.F
641!|| imp_fhht ../engine/source/implicit/imp_dyna.F
642!|| imp_fhht1 ../engine/source/implicit/imp_dyna.F
643!|| imp_setb ../engine/source/implicit/imp_setb.F
644!|| imp_setba ../engine/source/implicit/imp_setb.F
645!|| imp_setbp ../engine/source/implicit/imp_setb.F
646!|| rbe2_impr1 ../engine/source/constraints/general/rbe2/rbe2_imp0.F
647!|| rbe3_impr1 ../engine/source/constraints/general/rbe3/rbe3_imp0.F
648!|| rbe3_impr2 ../engine/source/constraints/general/rbe3/rbe3_imp0.F
649!|| rby_impr1 ../engine/source/constraints/general/rbody/rby_imp0.F
650!|| rby_impr2 ../engine/source/constraints/general/rbody/rby_imp0.F
651!|| spmd_sumf_v ../engine/source/mpi/implicit/imp_spmd.F
652!||--- uses -----------------------------------------------------
653!|| intbufdef_mod ../common_source/modules/interfaces/intbufdef_mod.f90
654!||====================================================================
655 SUBROUTINE upd_rhs(ICODT ,ICODR ,ISKEW ,IBFV ,XFRAME ,
656 1 RBY ,X ,SKEW ,LPBY ,NPBY ,
657 2 NRBYAC,IRBYAC,NINT2 ,IINT2 ,IPARI ,
658 3 INTBUF_TAB ,NDOF ,IDDL ,IKC ,
659 4 NDDL0 ,B ,IUPD ,INLOC ,LJ ,
660 5 A ,AR ,AC ,ACR ,NT_RW ,
661 6 IRFLAG,W_DDL ,NDDL ,R02 ,IDYNA ,
662 7 V ,VR ,MS ,IN ,IRBE3 ,
663 8 LRBE3 ,FRBE3 ,WEIGHT ,IRBE2 ,LRBE2 )
664C-----------------------------------------------
665C M o d u l e s
666C-----------------------------------------------
667 USE intbufdef_mod
668C-----------------------------------------------
669C I m p l i c i t T y p e s
670C-----------------------------------------------
671#include "implicit_f.inc"
672C-----------------------------------------------
673C C o m m o n B l o c k s
674C-----------------------------------------------
675#include "com01_c.inc"
676#include "com04_c.inc"
677#include "param_c.inc"
678C-----------------------------------------------
679C D u m m y A r g u m e n t s
680C-----------------------------------------------
681 INTEGER IBFV(NIFV,*),ICODT(*),ICODR(*),ISKEW(*),
682 . NINT2 ,IINT2(*),LJ(*),NDDL0,IUPD,
683 . INLOC(*),NT_RW,IRFLAG,W_DDL(*) ,NDDL,IDYNA
684 INTEGER LPBY(*),NPBY(NNPBY,*),NDOF(*),IDDL(*),IKC(*),
685 . IPARI(NPARI,*), NRBYAC,IRBYAC(*)
686 INTEGER WEIGHT(*),IRBE3(*),LRBE3(*),IRBE2(*),LRBE2(*)
687 my_real
688 . RBY(NRBY,*) ,X(3,*) ,SKEW(*),R02
689 my_real
690 . B(*) ,XFRAME(NXFRAME,*),A(3,*),AR(3,*),AC(3,*),ACR(3,*),
691 . v(3,*),vr(3,*),ms(*),in(*),frbe3(*)
692 TYPE(intbuf_struct_) INTBUF_TAB(*)
693C-----------------------------------------------
694C L o c a l V a r i a b l e s
695C-----------------------------------------------
696 INTEGER I,J,K,N,JI,JB,K1,IFLAG
697 my_real,
698 . DIMENSION(:,:),ALLOCATABLE :: MA,INA
699 my_real,
700 . DIMENSION(:),ALLOCATABLE :: B0
701C-------int2,RBE3,rby speciale (Fext seulement)----------
702 IF (IUPD==0) then
703 DO i=1,nint2
704 n=iint2(i)
705 CALL i2_impr1(ipari(1,n),intbuf_tab(n) ,
706 . x ,ndof ,iddl ,b )
707 ENDDO
708 IF (nrbe2>0) THEN
709 CALL rbe2_impr1(
710 1 irbe2 ,lrbe2 ,x ,skew ,ndof ,
711 2 iddl ,b ,weight)
712 ENDIF
713 IF (nrbe3>0) THEN
714 CALL rbe3_impr1(
715 1 irbe3 ,lrbe3 ,frbe3 ,x ,skew ,
716 2 ndof ,iddl ,b ,weight)
717 ENDIF
718 DO i=1,nrbyac
719 n=irbyac(i)
720 k1=irbyac(i+nrbykin)+1
721 CALL rby_impr1(x, rby(1,n),lpby(k1),npby(1,n),
722 1 ndof ,iddl ,b )
723 ENDDO
724 ENDIF
725C-------int2,rby speciale (elems deleted)----------
726 DO i=1,nint2
727 n=iint2(i)
728 CALL i2_impr2(ipari(1,n),intbuf_tab(n) ,ac ,acr ,
729 . x ,ndof ,iddl ,b )
730 ENDDO
731 IF (nrbe3>0) THEN
732 CALL rbe3_impr2(
733 1 irbe3 ,lrbe3 ,frbe3 ,x ,skew ,
734 2 ndof ,iddl ,b ,weight,ac ,
735 3 acr )
736 ENDIF
737 DO i=1,nrbyac
738 n=irbyac(i)
739 k1=irbyac(i+nrbykin)+1
740 CALL rby_impr2(x, rby(1,n),lpby(k1),npby(1,n),
741 1 ndof ,iddl ,b ,ac ,acr )
742 ENDDO
743 IF (idyna>0) THEN
744 ALLOCATE(ma(3,numnod))
745 IF (iroddl/=0) ALLOCATE(ina(3,numnod))
746 CALL imp_dynar(ma ,ina ,ms ,in ,a ,ar ,
747 . v ,vr )
748 DO i=1,nrbyac
749 n=irbyac(i)
750 CALL imp_dycrb(ina ,in ,vr ,npby(1,n),rby(1,n),
751 . weight,icodr ,iskew ,skew )
752 ENDDO
753 ENDIF
754 IF (irflag>0) THEN
755 CALL ext_rhs(icodt ,icodr ,iskew ,ibfv ,xframe ,
756 1 x ,skew ,ndof ,iddl ,ikc ,
757 2 nddl0 ,b ,inloc ,lj ,ac ,
758 3 acr ,nt_rw ,w_ddl ,nddl ,r02 )
759 ALLOCATE(b0(nddl0))
760 CALL imp_setb(a ,ar ,iddl ,ndof ,b0 )
761 CALL ext_rhs(icodt ,icodr ,iskew ,ibfv ,xframe ,
762 1 x ,skew ,ndof ,iddl ,ikc ,
763 2 nddl0 ,b0 ,inloc ,lj ,a ,
764 3 ar ,nt_rw ,w_ddl ,nddl ,r02 )
765 IF (idyna>0) THEN
766 CALL imp_setb(ma ,ina ,iddl ,ndof ,b0 )
767 CALL ext_rhs(icodt ,icodr ,iskew ,ibfv ,xframe ,
768 1 x ,skew ,ndof ,iddl ,ikc ,
769 2 nddl0 ,b0 ,inloc ,lj ,ma ,
770 3 ina ,nt_rw ,w_ddl ,nddl ,r02 )
771 END IF !(IDYNA>0) THEN
772 DEALLOCATE(b0)
773 ENDIF
774C -----add Fint---------------------
775 iflag = 1
776 CALL imp_setba(a ,ar ,iddl ,ndof ,b ,
777 1 iflag )
778 DO i = 1,numnod
779 a(1,i) = a(1,i) + ac(1,i)
780 a(2,i) = a(2,i) + ac(2,i)
781 a(3,i) = a(3,i) + ac(3,i)
782 ENDDO
783 IF (iroddl/=0) THEN
784 DO i = 1,numnod
785 ar(1,i) = ar(1,i) + acr(1,i)
786 ar(2,i) = ar(2,i) + acr(2,i)
787 ar(3,i) = ar(3,i) + acr(3,i)
788 ENDDO
789 ENDIF
790C+++
791 IF (idyna>0) THEN
792C-------------
793 CALL imp_fhht(nddl0 ,b )
794C -----add Ma----M,IN sont deja condenses-----------------
795 CALL imp_setba(ma ,ina ,iddl ,ndof ,b ,
796 1 iflag )
797 DO i = 1,numnod
798 a(1,i) = a(1,i) + ma(1,i)
799 a(2,i) = a(2,i) + ma(2,i)
800 a(3,i) = a(3,i) + ma(3,i)
801 ENDDO
802 IF (iroddl/=0) THEN
803 DO i = 1,numnod
804 ar(1,i) = ar(1,i) + ina(1,i)
805 ar(2,i) = ar(2,i) + ina(2,i)
806 ar(3,i) = ar(3,i) + ina(3,i)
807 ENDDO
808 ENDIF
809 DEALLOCATE(ma)
810 IF (iroddl/=0) DEALLOCATE(ina)
811 ENDIF
812C ---
813C
814C
815 CALL bc_impr1(icodt ,icodr ,iskew ,skew ,ndof ,
816 1 iddl ,b )
817C
818 IF(nfxvel>0) THEN
819 CALL fv_imprl(ibfv ,skew ,xframe ,lj ,iddl ,
820 1 ndof ,b )
821 ENDIF
822C
823 IF(nt_rw>0) CALL fv_rwlr0(iddl ,b )
824C
825 IF (nspmd>1) THEN
826C ----------B->dB---------------------
827 iflag = 0
828 CALL imp_setba(a ,ar ,iddl ,ndof ,b ,
829 1 iflag )
830 CALL condens_b(nddl0 ,ikc ,b )
831 CALL spmd_sumf_v(b )
832 CALL imp_setbp(a ,ar ,iddl ,ndof ,ikc ,
833 . inloc ,b )
834 ELSE
835 CALL condens_b(nddl0 ,ikc ,b )
836 ENDIF
837C-------------
838 IF (idyna>0) CALL imp_fhht1(nddl0 ,nddl ,b ,ikc )
839C
840 RETURN
841 END
842!||====================================================================
843!|| ext_rhs ../engine/source/implicit/upd_glob_k.F
844!||--- called by ------------------------------------------------------
845!|| upd_rhs ../engine/source/implicit/upd_glob_k.F
846!|| upd_rhs_fr ../engine/source/implicit/imp_solv.F
847!||--- calls -----------------------------------------------------
848!|| bc_impr1 ../engine/source/constraints/general/bcs/bc_imp0.F
849!|| condens_b ../engine/source/implicit/upd_glob_k.F
850!|| cp_real ../engine/source/implicit/produt_v.F
851!|| fv_imprl ../engine/source/constraints/general/impvel/fv_imp0.F
852!|| fv_rwlr0 ../engine/source/constraints/general/rwall/srw_imp.F
853!|| imp_setba ../engine/source/implicit/imp_setb.F
854!|| imp_setbp ../engine/source/implicit/imp_setb.F
855!|| produt_w ../engine/source/implicit/produt_v.F
856!|| spmd_sumf_v ../engine/source/mpi/implicit/imp_spmd.F
857!||====================================================================
858 SUBROUTINE ext_rhs(ICODT ,ICODR ,ISKEW ,IBFV ,XFRAME ,
859 1 X ,SKEW ,NDOF ,IDDL ,IKC ,
860 4 NDDL0 ,B ,INLOC ,LJ ,AC ,
861 5 ACR ,NT_RW ,W_DDL ,NDDL ,R02 )
862C-----------------------------------------------
863C I m p l i c i t T y p e s
864C-----------------------------------------------
865#include "implicit_f.inc"
866C-----------------------------------------------
867C C o m m o n B l o c k s
868C-----------------------------------------------
869#include "com01_c.inc"
870#include "com04_c.inc"
871#include "param_c.inc"
872C-----------------------------------------------
873C D u m m y A r g u m e n t s
874C-----------------------------------------------
875 INTEGER IBFV(NIFV,*),ICODT(*),ICODR(*),ISKEW(*),
876 . LJ(*),NDDL0,INLOC(*),NT_RW,W_DDL(*),NDDL
877 INTEGER NDOF(*),IDDL(*),IKC(*)
878 my_real
879 . x(3,*) ,skew(*),r02
880 my_real
881 . b(*) ,xframe(nxframe,*),ac(3,*),acr(3,*)
882C-----------------------------------------------
883C L o c a l V a r i a b l e s
884C-----------------------------------------------
885 INTEGER I,J,K,N,JI,JB,K1,IFLAG,ND
886 my_real
887 . LB(NDDL0) ,REXT
888C------------------------------------------
889 CALL CP_REAL(NDDL0,B,LB)
890 CALL BC_IMPR1(ICODT ,ICODR ,ISKEW ,SKEW ,NDOF ,
891 1 IDDL ,LB )
892 IF(NFXVEL>0) THEN
893 CALL FV_IMPRL(IBFV ,SKEW ,XFRAME ,LJ ,IDDL ,
894 1 NDOF ,LB )
895 ENDIF
896 IF(NT_RW>0) CALL FV_RWLR0(IDDL ,LB )
897C
898 IF (nspmd>1) THEN
899C ----------LB->dB---------------------
900 iflag = 0
901 CALL imp_setba(ac ,acr ,iddl ,ndof ,lb ,
902 1 iflag )
903 CALL condens_b(nddl0 ,ikc ,lb )
904 CALL spmd_sumf_v(lb )
905 CALL imp_setbp(ac ,acr ,iddl ,ndof ,ikc ,
906 . inloc ,lb )
907 ELSE
908 CALL condens_b(nddl0 ,ikc ,lb )
909 ENDIF
910C
911 CALL produt_w(nddl, lb, lb,w_ddl,rext )
912 r02=max(r02,rext )
913C
914 RETURN
915 END
916!||====================================================================
917!|| rer02 ../engine/source/implicit/upd_glob_k.F
918!||--- called by ------------------------------------------------------
919!|| imp_solv ../engine/source/implicit/imp_solv.F
920!||--- calls -----------------------------------------------------
921!|| imp_setbp ../engine/source/implicit/imp_setb.F
922!|| produt_w ../engine/source/implicit/produt_v.F
923!|| recukin ../engine/source/implicit/recudis.F
924!|| rer_int_v ../engine/source/implicit/upd_glob_k.F
925!|| spmd_sumf_v ../engine/source/mpi/implicit/imp_spmd.F
926!||--- uses -----------------------------------------------------
927!|| intbufdef_mod ../common_source/modules/interfaces/intbufdef_mod.F90
928!||====================================================================
929 SUBROUTINE rer02(RBY ,LPBY ,NPBY ,SKEW ,ISKEW ,
930 1 ITAB ,WEIGHT,MS ,IN ,
931 2 IBFV ,VEL ,ICODT,ICODR ,
932 3 NRBYAC,IRBYAC,NINT2 ,IINT2 ,IPARI ,
933 4 INTBUF_TAB ,NDOF ,D ,DR ,
934 5 X ,XFRAME,LJ ,IXR ,IXC ,
935 6 IXTG ,SH4TREE,SH3TREE,IRBE3 ,LRBE3,
936 7 FRBE3 ,IADK ,JDIK ,DIAG_K,LT_K ,
937 8 IDDL ,IKC ,INLOC ,NUM_IMP,NS_IMP,
938 9 NE_IMP,INDEX2,NDDL ,W_DDL ,A ,
939 A AR ,R02 ,IRBE2 ,LRBE2 ,X_C )
940C-----------------------------------------------
941C M o d u l e s
942C-----------------------------------------------
943 USE intbufdef_mod
944C-----------------------------------------------
945C I m p l i c i t T y p e s
946C-----------------------------------------------
947#include "implicit_f.inc"
948C-----------------------------------------------
949C C o m m o n B l o c k s
950C-----------------------------------------------
951#include "com01_c.inc"
952#include "com04_c.inc"
953#include "param_c.inc"
954#include "units_c.inc"
955#include "impl1_c.inc"
956#include "task_c.inc"
957C-----------------------------------------------
958C D u m m y A r g u m e n t s
959C-----------------------------------------------
960 INTEGER IBFV(NIFV,*),LJ(*)
961 INTEGER WEIGHT(*),LPBY(*),NPBY(NNPBY,*),ITAB(*),
962 . IPARI(NPARI,*), ISKEW(*),
963 . nrbyac,irbyac(*),nint2 ,iint2(*),ixr(*)
964 INTEGER NDOF(*),ICODT(*) ,ICODR(*),IXC(*),IXTG(*),
965 . SH4TREE(*), SH3TREE(*),IRBE3(*),LRBE3(*),
966 . IRBE2(*),LRBE2(*)
967 INTEGER NDDL,IADK(*),JDIK(*),IDDL(*),IKC(*) ,INLOC(*),
968 . NUM_IMP(*),NS_IMP(*),NE_IMP(*),INDEX2(*),W_DDL(*)
969 my_real
970 . RBY(NRBY,*) ,SKEW(*),IN(*),MS(*),
971 . VEL(LFXVELR,*), XFRAME(NXFRAME,*),FRBE3(*),X_C(*)
972 my_real
973 . X(3,*),D(3,*),DR(3,*),DIAG_K(*),LT_K(*),A(3,*),AR(3,*),R02
974
975 TYPE(INTBUF_STRUCT_) INTBUF_TAB(*)
976C----------------------------------------------
977C L o c a l V a r i a b l e s
978C-----------------------------------------------
979 INTEGER I,J,K,N,ND,NNZ
980 my_real
981 . U(NDDL),F(NDDL),R2
982C-------D=0 for free dof------------
983 DO N =1,numnod
984 a(1,n)= d(1,n)
985 a(2,n)= d(2,n)
986 a(3,n)= d(3,n)
987 ENDDO
988 IF (iroddl/=0) THEN
989 DO n =1,numnod
990 ar(1,n)= dr(1,n)
991 ar(2,n)= dr(2,n)
992 ar(3,n)= dr(3,n)
993 ENDDO
994 ENDIF
995 DO n =1,numnod
996 i=inloc(n)
997 DO j=1,min(3,ndof(i))
998 nd = iddl(i)+j
999 IF (ikc(nd)==0) a(j,i) =zero
1000 ENDDO
1001 DO j=4,ndof(i)
1002 nd = iddl(i)+j
1003 IF (ikc(nd)==0) ar(j-3,i) =zero
1004 ENDDO
1005 ENDDO
1006C
1007 IF(imp_lr > 0)THEN
1008 CALL recukin(rby ,lpby ,npby ,skew ,iskew ,
1009 1 itab ,weight,ms ,in ,
1010 2 ibfv ,vel ,icodt,icodr ,
1011 3 nrbyac,irbyac,nint2 ,iint2 ,ipari ,
1012 4 intbuf_tab,ndof ,a ,ar ,
1013 5 x_c ,xframe,lj ,ixr ,ixc ,
1014 6 ixtg ,sh4tree,sh3tree,irbe3 ,lrbe3,
1015 7 frbe3 ,irbe2 ,lrbe2 )
1016 ELSE
1017 CALL recukin(rby ,lpby ,npby ,skew ,iskew ,
1018 1 itab ,weight,ms ,in ,
1019 2 ibfv ,vel ,icodt,icodr ,
1020 3 nrbyac,irbyac,nint2 ,iint2 ,ipari ,
1021 4 intbuf_tab ,ndof ,a ,ar ,
1022 5 x ,xframe,lj ,ixr ,ixc ,
1023 6 ixtg ,sh4tree,sh3tree,irbe3 ,lrbe3,
1024 7 frbe3 ,irbe2 ,lrbe2 )
1025 END IF !(IMP_LR > 0)THEN
1026C
1027 CALL rer_int_v(
1028 1 ipari ,intbuf_tab,num_imp ,ns_imp ,ne_imp ,
1029 2 x ,a )
1030 DO i =1,nddl
1031 u(i) =zero
1032 ENDDO
1033 CALL imp_setbp(a ,ar ,iddl ,ndof ,ikc ,
1034 . inloc ,u )
1035 DO i =1,nddl
1036 f(i)=diag_k(i)*u(i)
1037 ENDDO
1038 IF (nspmd > 1)CALL spmd_sumf_v(f )
1039c NNZ = IADK(NDDL+1)-IADK(1)
1040c CALL MAV_LT(
1041c 1 NDDL ,NNZ ,IADK ,JDIK ,DIAG_K ,
1042c 2 LT_K ,U ,F )
1043 CALL produt_w( nddl ,f ,f ,w_ddl , r2)
1044 IF (r2>r02.AND.imconv>=0.AND.
1045 . (lprint/=0.OR.nprint/=0)) THEN
1046 IF(ispmd==0) THEN
1047 WRITE(iout,1006)sqrt(r2)
1048 IF(nprint<0) THEN
1049 WRITE(istdo,1006)sqrt(r2)
1050 ENDIF
1051 ENDIF
1052 ENDIF
1053 r02= max(r2,r02)
1054C
1055 1006 FORMAT(3x,'DUE TO THE CONTACT,',
1056 . ' REFERENCE RESIDUAL(F) NORM HAS BEEN CHANGED TO:',e11.4)
1057 RETURN
1058 END
1059!||====================================================================
1060!|| rer_int_v ../engine/source/implicit/upd_glob_k.F
1061!||--- called by ------------------------------------------------------
1062!|| rer02 ../engine/source/implicit/upd_glob_k.f
1063!||--- calls -----------------------------------------------------
1064!|| ud_int11 ../engine/source/implicit/upd_glob_k.F
1065!|| ud_int5 ../engine/source/implicit/upd_glob_k.F
1066!|| ud_int7 ../engine/source/implicit/upd_glob_k.F
1067!||--- uses -----------------------------------------------------
1068!|| intbufdef_mod ../common_source/modules/interfaces/intbufdef_mod.F90
1069!||====================================================================
1070 SUBROUTINE rer_int_v(
1071 1 IPARI ,INTBUF_TAB,NUM_IMP ,NS_IMP ,NE_IMP ,
1072 2 X ,UD )
1073C-----------------------------------------------
1074C M o d u l e s
1075C-----------------------------------------------
1076 USE intbufdef_mod
1077C-----------------------------------------------
1078C I m p l i c i t T y p e s
1079C-----------------------------------------------
1080#include "implicit_f.inc"
1081C-----------------------------------------------
1082C C o m m o n B l o c k s
1083C-----------------------------------------------
1084#include "com04_c.inc"
1085#include "param_c.inc"
1086C-----------------------------------------------
1087C D u m m y A r g u m e n t s
1088C-----------------------------------------------
1089 INTEGER IPARI(NPARI,*),NUM_IMP(*),
1090 . NS_IMP(*),NE_IMP(*)
1091C REAL
1092 my_real
1093 . X(3,*),UD(3,*)
1094 TYPE(intbuf_struct_) INTBUF_TAB(*)
1095C-----------------------------------------------
1096C L o c a l V a r i a b l e s
1097C-----------------------------------------------
1098 INTEGER NIN,NTY,NSN
1099 INTEGER I,J,K,L,NDOFI,N,IAD,
1100 . idrec,insv,ins11,nrts
1101C-----------------------------------------------
1102 iad=1
1103 DO nin=1,ninter
1104 nty =ipari(7,nin)
1105 nsn =ipari(5,nin)
1106C
1107 IF(nty==3)THEN
1108 ELSEIF(nty==4)THEN
1109 ELSEIF(nty==5)THEN
1110 CALL ud_int5(
1111 1 ipari(1,nin),intbuf_tab(nin) ,x ,
1112 2 num_imp(nin),ns_imp(iad),ne_imp(iad) ,ud )
1113 iad=iad+num_imp(nin)
1114 ENDIF
1115 ENDDO
1116C
1117 DO nin=1,ninter
1118 nty =ipari(7,nin)
1119 nsn =ipari(5,nin)
1120C
1121 IF(nty==3)THEN
1122 ELSEIF(nty==4)THEN
1123 ELSEIF(nty==5)THEN
1124 ELSEIF(nty==6)THEN
1125C------------verify effect, if R02 is too high-----
1126 ELSEIF(nty==7.OR.nty==10.OR.nty==24)THEN
1127C
1128 CALL ud_int7(
1129 1 num_imp(nin),ns_imp(iad),ne_imp(iad),intbuf_tab(nin)%IRECTM,
1130 . intbuf_tab(nin)%NSV,
1131 2 nsn ,x ,ud )
1132 iad=iad+num_imp(nin)
1133 ELSEIF(nty==11)THEN
1134C
1135 nrts =ipari(3,nin)
1136 CALL ud_int11(
1137 1 num_imp(nin),ns_imp(iad),ne_imp(iad),intbuf_tab(nin)%IRECTS,
1138 2 intbuf_tab(nin)%IRECTM,nsn ,x ,ud )
1139 iad=iad+num_imp(nin)
1140 ENDIF
1141 ENDDO
1142C----6---------------------------------------------------------------7---------8
1143 RETURN
1144 END
1145C --------- ----
1146!||====================================================================
1147!|| ud_int7 ../engine/source/implicit/upd_glob_k.F
1148!||--- called by ------------------------------------------------------
1149!|| rer_int_v ../engine/source/implicit/upd_glob_k.F
1150!||--- calls -----------------------------------------------------
1151!|| inihi7 ../engine/source/implicit/upd_glob_k.F
1152!||====================================================================
1153 SUBROUTINE ud_int7(
1154 1 JLT ,NS_IMP ,NE_IMP ,IRECT ,NSV ,
1155 2 NSN ,X ,UD )
1156C-----------------------------------------------
1157C I m p l i c i t T y p e s
1158C-----------------------------------------------
1159#include "implicit_f.inc"
1160C-----------------------------------------------
1161C D u m m y A r g u m e n t s
1162C-----------------------------------------------
1163 INTEGER JLT,NS_IMP(*),NE_IMP(*),IRECT(4,*),NSV(*),NSN
1164C REAL
1165 my_real
1166 . x(3,*),ud(3,*)
1167C-----------------------------------------------
1168C L o c a l V a r i a b l e s
1169C-----------------------------------------------
1170 INTEGER I,J,N,N1,N2,NE,IG,M(4)
1171 my_real
1172 . h(4),us(3)
1173C-------Missed SPMD----------------------------------------
1174 DO i = 1, jlt
1175C--------secnd node-----add case Ud_s >0
1176 ig = ns_imp(i)
1177 IF (ig<=nsn) THEN
1178 n1 = nsv(ig)
1179 ne=ne_imp(i)
1180 DO j = 1, 4
1181 m(j)= irect(j,ne)
1182 ENDDO
1183 CALL inihi7(n1,m,x ,h)
1184 IF ((abs(ud(1,n1))+abs(ud(2,n1))+abs(ud(3,n1)))>zero) THEN
1185 DO j = 1, 4
1186 us(1) = h(j)*ud(1,n1)
1187 us(2) = h(j)*ud(2,n1)
1188 us(3) = h(j)*ud(3,n1)
1189 IF (abs(us(1))>abs(ud(1,m(j)))) ud(1,m(j))= us(1)
1190 IF (abs(us(2))>abs(ud(2,m(j)))) ud(2,m(j))= us(2)
1191 IF (abs(us(3))>abs(ud(3,m(j)))) ud(3,m(j))= us(3)
1192 ENDDO
1193 ELSE
1194 us(1) = zero
1195 us(2) = zero
1196 us(3) = zero
1197 DO j = 1, 4
1198 us(1) = us(1)+ h(j)*ud(1,m(j))
1199 us(2) = us(2)+ h(j)*ud(2,m(j))
1200 us(3) = us(3)+ h(j)*ud(3,m(j))
1201 ENDDO
1202 IF (abs(us(1))>abs(ud(1,n1))) ud(1,n1)= us(1)
1203 IF (abs(us(2))>abs(ud(2,n1))) ud(2,n1)= us(2)
1204 IF (abs(us(3))>abs(ud(3,n1))) ud(3,n1)= us(3)
1205 END IF
1206 ENDIF
1207 ENDDO
1208C----6---------------------------------------------------------------7---------8
1209 RETURN
1210 END
1211!||====================================================================
1212!|| ud_int11 ../engine/source/implicit/upd_glob_k.F
1213!||--- called by ------------------------------------------------------
1214!|| rer_int_v ../engine/source/implicit/upd_glob_k.F
1215!||--- calls -----------------------------------------------------
1216!|| inihi11 ../engine/source/implicit/upd_glob_k.F
1217!||====================================================================
1218 SUBROUTINE ud_int11(
1219 1 JLT ,NS_IMP ,NE_IMP ,IRECTS ,IRECTM ,
1220 2 NSN ,X ,UD )
1221C-----------------------------------------------
1222C I m p l i c i t T y p e s
1223C-----------------------------------------------
1224#include "implicit_f.inc"
1225C-----------------------------------------------
1226C D u m m y A r g u m e n t s
1227C-----------------------------------------------
1228 INTEGER JLT,NS_IMP(*),NE_IMP(*),IRECTS(2,*),IRECTM(2,*),NSN
1229C REAL
1230 my_real
1231 . x(3,*),ud(3,*)
1232C-----------------------------------------------
1233C L o c a l V a r i a b l e s
1234C-----------------------------------------------
1235 INTEGER I,J,N,N1,N2,NE,IG,M1,M2,NS(2),NM(2)
1236 my_real
1237 . hs(2) ,hm(2),us(3,2),vm
1238C-----------------------------------------------
1239 DO i = 1, jlt
1240C--------secnd node-----
1241 ig = ns_imp(i)
1242 IF (ig<=nsn) THEN
1243 ns(1) = irects(1,ig)
1244 ns(2) = irects(2,ig)
1245 ne=ne_imp(i)
1246 nm(1) = irectm(1,ne)
1247 nm(2) = irectm(2,ne)
1248 CALL inihi11(ns, nm, x ,hs ,hm)
1249 DO j = 1, 3
1250 vm = hm(1)*ud(j,nm(1))+hm(2)*ud(j,nm(2))
1251 us(j,1)=hs(1)*vm
1252 us(j,2)=hs(2)*vm
1253 IF (abs(us(j,1))>abs(ud(j,ns(1)))) ud(j,ns(1))= us(j,1)
1254 IF (abs(us(j,2))>abs(ud(j,ns(2)))) ud(j,ns(2))= us(j,2)
1255 ENDDO
1256 ENDIF
1257 ENDDO
1258C----6---------------------------------------------------------------7---------8
1259 RETURN
1260 END
1261!||====================================================================
1262!|| inihi7 ../engine/source/implicit/upd_glob_k.F
1263!||--- called by ------------------------------------------------------
1264!|| ud_int7 ../engine/source/implicit/upd_glob_k.F
1265!||====================================================================
1266 SUBROUTINE inihi7(NS,IRECT,X ,H)
1267C
1268C-----------------------------------------------
1269C I m p l i c i t T y p e s
1270C-----------------------------------------------
1271#include "implicit_f.inc"
1272C-----------------------------------------------
1273C D u m m y A r g u m e n t s
1274C-----------------------------------------------
1275 INTEGER NS,IRECT(4)
1276C REAL
1277 my_real
1278 . X(3,*) ,H(4)
1279C-----------------------------------------------
1280C L o c a l V a r i a b l e s
1281C-----------------------------------------------
1282C REAL
1283 INTEGER N1,N2,N3,N4
1284 my_real
1285 . x0, y0, z0, xl1, xl2, xl3, xl4, yy1, yy2, yy3, yy4,
1286 . zz1, zz2, zz3, zz4, xi1, xi2, xi3, xi4, yi1, yi2, yi3, yi4,
1287 . zi1, zi2, zi3, zi4, xn1, yn1, zn1, xn2, yn2, zn2, xn3, yn3,
1288 . zn3, xn4, yn4, zn4, an, area, a12, a23, a34, a41, b12, b23,
1289 . b34, b41, ab1, ab2, tp, tm, sp, sm, x1,x2,x3,x4,
1290 . y1,y2,y3,y4,z1,z2,z3,z4,xi,yi,zi,nx,ny,nz,alp,ssc,ttc
1291C-----------------------------------------------
1292 alp = twoem2
1293 h(1) = zero
1294 h(2) = zero
1295 h(3) = zero
1296 h(4) = zero
1297 n1 = irect(1)
1298 n2 = irect(2)
1299 n3 = irect(3)
1300 n4 = irect(4)
1301 x1=x(1,n1)
1302 x2=x(1,n2)
1303 x3=x(1,n3)
1304 x4=x(1,n4)
1305 y1=x(2,n1)
1306 y2=x(2,n2)
1307 y3=x(2,n3)
1308 y4=x(2,n4)
1309 z1=x(3,n1)
1310 z2=x(3,n2)
1311 z3=x(3,n3)
1312 z4=x(3,n4)
1313 xi=x(1,ns)
1314 yi=x(2,ns)
1315 zi=x(3,ns)
1316C
1317 x0 = fourth*(x1+x2+x3+x4)
1318 y0 = fourth*(y1+y2+y3+y4)
1319 z0 = fourth*(z1+z2+z3+z4)
1320C
1321 xl1 = x1-x0
1322 xl2 = x2-x0
1323 xl3 = x3-x0
1324 xl4 = x4-x0
1325 yy1 = y1-y0
1326 yy2 = y2-y0
1327 yy3 = y3-y0
1328 yy4 = y4-y0
1329 zz1 = z1-z0
1330 zz2 = z2-z0
1331 zz3 = z3-z0
1332 zz4 = z4-z0
1333C
1334 xi1 = x1-xi
1335 xi2 = x2-xi
1336 xi3 = x3-xi
1337 xi4 = x4-xi
1338 yi1 = y1-yi
1339 yi2 = y2-yi
1340 yi3 = y3-yi
1341 yi4 = y4-yi
1342 zi1 = z1-zi
1343 zi2 = z2-zi
1344 zi3 = z3-zi
1345 zi4 = z4-zi
1346C
1347 xn1 = yy1*zz2 - yy2*zz1
1348 yn1 = zz1*xl2 - zz2*xl1
1349 zn1 = xl1*yy2 - xl2*yy1
1350 nx=xn1
1351 ny=yn1
1352 nz=zn1
1353C
1354 xn2 = yy2*zz3 - yy3*zz2
1355 yn2 = zz2*xl3 - zz3*xl2
1356 zn2 = xl2*yy3 - xl3*yy2
1357 nx=nx+xn2
1358 ny=ny+yn2
1359 nz=nz+zn2
1360C
1361 xn3 = yy3*zz4 - yy4*zz3
1362 yn3 = zz3*xl4 - zz4*xl3
1363 zn3 = xl3*yy4 - xl4*yy3
1364 nx=nx+xn3
1365 ny=ny+yn3
1366 nz=nz+zn3
1367C
1368 xn4 = yy4*zz1 - yy1*zz4
1369 yn4 = zz4*xl1 - zz1*xl4
1370 zn4 = xl4*yy1 - xl1*yy4
1371 nx=nx+xn3
1372 ny=ny+yn3
1373 nz=nz+zn3
1374C
1375 an= max(em20,sqrt(nx*nx+ny*ny+nz*nz))
1376 IF(an<=em19) RETURN
1377 nx=nx/an
1378 ny=ny/an
1379 nz=nz/an
1380 area=half*an
1381C
1382 a12=(nx*xn1+ny*yn1+nz*zn1)
1383 a23=(nx*xn2+ny*yn2+nz*zn2)
1384 a34=(nx*xn3+ny*yn3+nz*zn3)
1385 a41=(nx*xn4+ny*yn4+nz*zn4)
1386C
1387 xn1 = yi1*zi2 - yi2*zi1
1388 yn1 = zi1*xi2 - zi2*xi1
1389 zn1 = xi1*yi2 - xi2*yi1
1390 b12=(nx*xn1+ny*yn1+nz*zn1)
1391C
1392 xn2 = yi2*zi3 - yi3*zi2
1393 yn2 = zi2*xi3 - zi3*xi2
1394 zn2 = xi2*yi3 - xi3*yi2
1395 b23=(nx*xn2+ny*yn2+nz*zn2)
1396C
1397 xn3 = yi3*zi4 - yi4*zi3
1398 yn3 = zi3*xi4 - zi4*xi3
1399 zn3 = xi3*yi4 - xi4*yi3
1400 b34=(nx*xn3+ny*yn3+nz*zn3)
1401C
1402 xn4 = yi4*zi1 - yi1*zi4
1403 yn4 = zi4*xi1 - zi1*xi4
1404 zn4 = xi4*yi1 - xi1*yi4
1405 b41=(nx*xn4+ny*yn4+nz*zn4)
1406C
1407 ab1=a23*b41
1408 ab2=b23*a41
1409C
1410 IF(abs(ab1+ab2)/area>em10)THEN
1411 ssc=(ab1-ab2)/(ab1+ab2)
1412 ELSE
1413 ssc=zero
1414 ENDIF
1415 IF(abs(a34/area)>em10)THEN
1416 ab1=b12*a34
1417 ab2=b34*a12
1418 ttc=(ab1-ab2)/(ab1+ab2)
1419 ELSE
1420 ttc=(b12-a12)/a12
1421 IF(b23<=zero.AND.b41<=zero)THEN
1422 IF(-b23/a12<=alp.AND.-b41/a12<=alp)ssc=zero
1423 ELSEIF(b23<=zero)THEN
1424 IF(-b23/a12<=alp)ssc=one
1425 ELSEIF(b41==zero)THEN
1426 IF(b41<zero.AND.-b41/a12<=alp)ssc=-one
1427 ENDIF
1428 ENDIF
1429C
1430 IF(abs(ssc)>one+alp.OR.abs(ttc)>one+alp) RETURN
1431 IF(abs(ssc)>one)ssc=ssc/abs(ssc)
1432 IF(abs(ttc)>one)ttc=ttc/abs(ttc)
1433C
1434 tp=fourth*(one+ttc)
1435 tm=fourth*(one-ttc)
1436C
1437 sp=one+ssc
1438 sm=one-ssc
1439 h(1)=tm*sm
1440 h(2)=tm*sp
1441 h(3)=tp*sp
1442 h(4)=tp*sm
1443C
1444 RETURN
1445 END
1446!||====================================================================
1447!|| inihi11 ../engine/source/implicit/upd_glob_k.f
1448!||--- called by ------------------------------------------------------
1449!|| ud_int11 ../engine/source/implicit/upd_glob_k.F
1450!||====================================================================
1451 SUBROUTINE inihi11(NS, NM, X ,HS ,HM)
1452C
1453C-----------------------------------------------
1454C I m p l i c i t T y p e s
1455C-----------------------------------------------
1456#include "implicit_f.inc"
1457C-----------------------------------------------
1458C D u m m y A r g u m e n t s
1459C-----------------------------------------------
1460 INTEGER NS(2),NM(2)
1461C REAL
1462 my_real
1463 . x(3,*) ,hs(2),hm(2)
1464C-----------------------------------------------
1465C L o c a l V a r i a b l e s
1466C-----------------------------------------------
1467C REAL
1468 my_real
1469 . xs12,ys12,zs12,xm12,ym12,zm12,xa,xb,
1470 . xs2,xm2,xsm,xs2m2,ys2,ym2,ysm,ys2m2,zs2,zm2,zsm,zs2m2,
1471 . xx,yy,zz,als,alm,det
1472C-----------------------------------------------
1473 hs(1) = zero
1474 hs(2) = zero
1475 hm(1) = zero
1476 hm(2) = zero
1477 xs12 = x(1,ns(2))-x(1,ns(1))
1478 ys12 = x(2,ns(2))-x(2,ns(1))
1479 zs12 = x(3,ns(2))-x(3,ns(1))
1480 xs2 = xs12*xs12 + ys12*ys12 + zs12*zs12
1481 xm12 = x(1,nm(2))-x(1,nm(1))
1482 ym12 = x(2,nm(2))-x(2,nm(1))
1483 zm12 = x(3,nm(2))-x(3,nm(1))
1484 xm2 = xm12*xm12 + ym12*ym12 + zm12*zm12
1485 xsm = - (xs12*xm12 + ys12*ym12 + zs12*zm12)
1486 xs2m2 = x(1,nm(2))-x(1,ns(2))
1487 ys2m2 = x(2,nm(2))-x(2,ns(2))
1488 zs2m2 = x(3,nm(2))-x(3,ns(2))
1489C
1490 xa = xs12*xs2m2 + ys12*ys2m2 + zs12*zs2m2
1491 xb = -xm12*xs2m2 - ym12*ys2m2 - zm12*zs2m2
1492 det = xm2*xs2 - xsm*xsm
1493 det = max(em20,det)
1494C
1495 hs(1) = (xb*xsm-xa*xm2) / det
1496 hm(1) = (xa*xsm-xb*xs2) / det
1497C
1498 xs2 = max(xs2,em20)
1499 xm2 = max(xm2,em20)
1500 IF(hm(1)<zero)THEN
1501 hm(1) = zero
1502 hs(1) = -xa / xs2
1503 ELSEIF(hm(1)>one)THEN
1504 hm(1) = one
1505 hs(1) = -(xa + xsm) / xs2
1506 ENDIF
1507C
1508 IF(hs(1)<zero)THEN
1509 hs(1) = zero
1510 hm(1) = -xb / xm2
1511 ELSEIF(hs(1)>one)THEN
1512 hs(1) = one
1513 hm(1) = -(xb + xsm) / xm2
1514 ENDIF
1515C
1516 hm(1) = min(one,hm(1))
1517 hm(1) = max(zero,hm(1))
1518C
1519 hs(2) = one -hs(1)
1520 hm(2) = one -hm(1)
1521C
1522 RETURN
1523 END
1524C --------- ----
1525!||====================================================================
1526!|| ud_int5 ../engine/source/implicit/upd_glob_k.F
1527!||--- called by ------------------------------------------------------
1528!|| rer_int_v ../engine/source/implicit/upd_glob_k.F
1529!||--- calls -----------------------------------------------------
1530!|| i3cst3 ../engine/source/interfaces/inter3d/i3cst3.F
1531!|| i3dis3 ../engine/source/interfaces/inter3d/i3dis3.F
1532!|| i3gap3 ../engine/source/interfaces/inter3d/i3gap3.F
1533!|| i5cork3 ../engine/source/interfaces/inter3d/i5cork3.F
1534!|| ud_intg5 ../engine/source/implicit/upd_glob_k.F
1535!||--- uses -----------------------------------------------------
1536!|| intbufdef_mod ../common_source/modules/interfaces/intbufdef_mod.F90
1537!||====================================================================
1538 SUBROUTINE ud_int5(
1539 1 IPARI ,INTBUF_TAB,X ,NUM_IMP,
1540 2 CAND_N,CAND_E ,UD )
1541C-----------------------------------------------
1542C M o d u l e s
1543C-----------------------------------------------
1544 USE intbufdef_mod
1545C----6---------------------------------------------------------------7---------8
1546C I m p l i c i t T y p e s
1547C-----------------------------------------------
1548#include "implicit_f.inc"
1549C-----------------------------------------------
1550C G l o b a l P a r a m e t e r s
1551C-----------------------------------------------
1552#include "mvsiz_p.inc"
1553C-----------------------------------------------
1554C C o m m o n B l o c k s
1555C-----------------------------------------------
1556#include "param_c.inc"
1557#include "vect01_c.inc"
1558C-----------------------------------------------------------------
1559C D u m m y A r g u m e n t s
1560C-----------------------------------------------
1561 INTEGER IPARI(*)
1562 INTEGER NUM_IMP,CAND_N(*),CAND_E(*)
1563C REAL
1564 my_real
1565 . x(3,*),ud(3,*)
1566
1567 TYPE(intbuf_struct_) INTBUF_TAB
1568C-----------------------------------------------
1569C L o c a l V a r i a b l e s
1570C-----------------------------------------------
1571 INTEGER I,IGAP, INACTI, IFQ, MFROT, IGSTI,INTY
1572 INTEGER JX1(MVSIZ), JX2(MVSIZ), JX3(MVSIZ), JX4(MVSIZ),
1573 . nsvg(mvsiz), i3n ,igimp
1574 INTEGER, DIMENSION(MVSIZ) :: IX1,IX2,IX3,IX4
1575C REAL
1576 my_real
1577 . startt, fric, gap, stopt,dist(mvsiz)
1578 my_real, DIMENSION(MVSIZ) :: x1,x2,x3,x4,xi
1579 my_real, DIMENSION(MVSIZ) :: y1,y2,y3,y4,yi
1580 my_real, DIMENSION(MVSIZ) :: z1,z2,z3,z4,zi
1581 my_real, DIMENSION(MVSIZ) :: xface,n1,n2,n3
1582 my_real, DIMENSION(MVSIZ) :: ssc,ttc,area,thk,alp
1583 my_real, DIMENSION(MVSIZ) :: x0,y0,z0,ans
1584 my_real, DIMENSION(MVSIZ) :: xx1,xx2,xx3,xx4
1585 my_real, DIMENSION(MVSIZ) :: yy1,yy2,yy3,yy4
1586 my_real, DIMENSION(MVSIZ) :: zz1,zz2,zz3,zz4
1587 my_real, DIMENSION(MVSIZ) :: xi1,xi2,xi3,xi4
1588 my_real, DIMENSION(MVSIZ) :: yi1,yi2,yi3,yi4
1589 my_real, DIMENSION(MVSIZ) :: zi1,zi2,zi3,zi4
1590 my_real, DIMENSION(MVSIZ) :: xn1,xn2,xn3,xn4
1591 my_real, DIMENSION(MVSIZ) :: yn1,yn2,yn3,yn4
1592 my_real, DIMENSION(MVSIZ) :: zn1,zn2,zn3,zn4
1593 my_real, DIMENSION(MVSIZ) :: xp,yp,zp
1594 my_real, DIMENSION(MVSIZ) :: h1,h2,h3,h4
1595C-----------------------------------------------
1596 inty = ipari(7)
1597 gap =intbuf_tab%VARIABLES(2)
1598 DO nft = 0 , num_imp - 1 , nvsiz
1599 lft=1
1600 llt = min( nvsiz, num_imp - nft )
1601 CALL i5cork3(
1602 1 x, intbuf_tab%IRECTM,intbuf_tab%MSR, intbuf_tab%NSV,
1603 2 intbuf_tab%IRTLM, cand_n(nft+1), cand_e(nft+1), nsvg,
1604 3 jx1, jx2, jx3, jx4,
1605 4 x1, x2, x3, x4,
1606 5 y1, y2, y3, y4,
1607 6 z1, z2, z3, z4,
1608 7 xface, xi, yi, zi,
1609 8 ix1, ix2, ix3, ix4,
1610 9 lft, llt, nft)
1611 CALL i3cst3(
1612 1 x1, x2, x3, x4,
1613 2 xi, y1, y2, y3,
1614 3 y4, yi, z1, z2,
1615 4 z3, z4, zi, xface,
1616 5 n1, n2, n3, ssc,
1617 6 ttc, x0, y0, z0,
1618 7 xx1, xx2, xx3, xx4,
1619 8 yy1, yy2, yy3, yy4,
1620 9 zz1, zz2, zz3, zz4,
1621 a xi1, xi2, xi3, xi4,
1622 b yi1, yi2, yi3, yi4,
1623 c zi1, zi2, zi3, zi4,
1624 d xn1, xn2, xn3, xn4,
1625 e yn1, yn2, yn3, yn4,
1626 f zn1, zn2, zn3, zn4,
1627 g area, lft, llt)
1628 CALL i3gap3(
1629 1 gap, area, thk, alp,
1630 2 lft, llt)
1631 CALL i3dis3(
1632 1 igimp, inty, dist, x1,
1633 2 x2, x3, x4, xi,
1634 3 y1, y2, y3, y4,
1635 4 yi, z1, z2, z3,
1636 5 z4, zi, xface, n1,
1637 6 n2, n3, ssc, ttc,
1638 7 alp, ans, xp, yp,
1639 8 zp, h1, h2, h3,
1640 9 h4, lft, llt)
1641 CALL ud_intg5(lft ,llt ,nsvg ,jx1 ,jx2 ,
1642 . jx3 ,jx4 ,ud ,h1 ,h2 ,
1643 2 h3 ,h4 ,xface)
1644 END DO
1645C----6---------------------------------------------------------------7---------8
1646 RETURN
1647 END
1648C --------- ----
1649!||====================================================================
1650!|| ud_intg5 ../engine/source/implicit/upd_glob_k.F
1651!||--- called by ------------------------------------------------------
1652!|| ud_int5 ../engine/source/implicit/upd_glob_k.F
1653!||====================================================================
1654 SUBROUTINE ud_intg5(LFT ,LLT ,NSVG ,JX1 ,JX2 ,
1655 . JX3 ,JX4 ,UD ,H1 ,H2 ,
1656 2 H3 ,H4 ,XFACE)
1657C----6---------------------------------------------------------------7---------8
1658C I m p l i c i t T y p e s
1659C-----------------------------------------------
1660#include "implicit_f.inc"
1661C-----------------------------------------------
1662C G l o b a l P a r a m e t e r s
1663C-----------------------------------------------
1664#include "mvsiz_p.inc"
1665C-----------------------------------------------
1666C D u m m y A r g u m e n t s
1667C-----------------------------------------------
1668 INTEGER LFT ,LLT ,NSVG(*) ,JX1(*), JX2(*),
1669 . JX3(*), JX4(*)
1670C REAL
1671 my_real
1672 . ud(3,*)
1673 my_real, DIMENSION(MVSIZ), INTENT(IN) :: h1,h2,h3,h4,xface
1674C-----------------------------------------------
1675C L o c a l V a r i a b l e s
1676C-----------------------------------------------
1677 INTEGER I, IL, IG, L, NN,M(4),J,NS
1678 my_real
1679 . H(4),US(3)
1680C-----------------------------------------------
1681 DO I = lft ,llt
1682C--------secnd node-----
1683 ns = nsvg(i)
1684 m(1)= jx1(i)
1685 m(2)= jx2(i)
1686 m(3)= jx3(i)
1687 m(4)= jx4(i)
1688 h(1)= h1(i)*xface(i)
1689 h(2)= h2(i)*xface(i)
1690 h(3)= h3(i)*xface(i)
1691 h(4)= h4(i)*xface(i)
1692 us(1) = zero
1693 us(2) = zero
1694 us(3) = zero
1695 DO j = 1, 4
1696 us(1) = us(1)+ h(j)*ud(1,m(j))
1697 us(2) = us(2)+ h(j)*ud(2,m(j))
1698 us(3) = us(3)+ h(j)*ud(3,m(j))
1699 ENDDO
1700 IF (abs(us(1))>abs(ud(1,ns))) ud(1,ns)= us(1)
1701 IF (abs(us(2))>abs(ud(2,ns))) ud(2,ns)= us(2)
1702 IF (abs(us(3))>abs(ud(3,ns))) ud(3,ns)= us(3)
1703 ENDDO
1704C----6---------------------------------------------------------------7---------8
1705 RETURN
1706 END
1707!||====================================================================
1708!|| ini_dofspc ../engine/source/implicit/upd_glob_k.F
1709!||--- called by ------------------------------------------------------
1710!|| upd_glob_k ../engine/source/implicit/upd_glob_k.F
1711!||--- calls -----------------------------------------------------
1712!|| colinear3 ../engine/source/implicit/upd_glob_k.F
1713!||--- uses -----------------------------------------------------
1714!|| intbufdef_mod ../common_source/modules/interfaces/intbufdef_mod.F90
1715!||====================================================================
1716 SUBROUTINE ini_dofspc(
1717 1 NPBY ,LPBY ,NRBYAC ,IRBYAC ,NINT2 ,
1718 2 IINT2 ,IPARI ,INTBUF_TAB,NDOF ,IRBE3 ,
1719 3 LRBE3 ,IRBE2 ,LRBE2 ,X ,DOFSPC )
1720C-----------------------------------------------
1721C M o d u l e s
1722C-----------------------------------------------
1723 USE intbufdef_mod
1724C-----------------------------------------------
1725C I m p l i c i t T y p e s
1726C-----------------------------------------------
1727#include "implicit_f.inc"
1728C-----------------------------------------------
1729C C o m m o n B l o c k s
1730C-----------------------------------------------
1731#include "param_c.inc"
1732#include "com04_c.inc"
1733C-----------------------------------------------
1734C D u m m y A r g u m e n t s
1735C-----------------------------------------------
1736 INTEGER NPBY(NNPBY,*),NRBYAC,IRBYAC(*),LPBY(*),DOFSPC(*)
1737 INTEGER NINT2,IINT2(*),IPARI(NPARI,*),NDOF(*),
1738 . IRBE3(NRBE3L,*) ,LRBE3(*),IRBE2(NRBE2L,*),LRBE2(*)
1739C REAL
1740 my_real
1741 . x(3,*)
1742
1743 TYPE(intbuf_struct_) INTBUF_TAB(*)
1744C-----------------------------------------------
1745C L o c a l V a r i a b l e s
1746C-----------------------------------------------
1747 INTEGER NMN,JI,K10,K11,K12,K13,K14,J,K
1748 INTEGER I,N,M,L,NS,ID,NM,NSN,IAD,ICOL,N1,N2,N3
1749C-----------------------------------------------
1750 DO i=1,numnod
1751 dofspc(i)=ndof(i)
1752 END DO
1753c RETURN
1754C
1755 DO i=1,nrbyac
1756 n=irbyac(i)
1757 k=irbyac(i+nrbykin)
1758 m=npby(1,n)
1759 IF (m<=0) cycle
1760 dofspc(m)=0
1761 nsn =npby(2,n)
1762 icol = 1
1763 DO n =1,nsn
1764 ns=lpby(k+n)
1765 IF (ndof(ns)>3) THEN
1766 icol = 0
1767 cycle
1768 ENDIF
1769 END DO
1770 IF (icol==1) THEN
1771 IF (nsn<=2) THEN
1772 dofspc(m)=ndof(m)
1773 ELSE
1774 CALL colinear3(lpby(k+1),lpby(k+2),lpby(k+3),x,icol)
1775 IF (icol==1) THEN
1776 n1=lpby(k+2)
1777 n2=lpby(k+3)
1778 DO n =4,nsn
1779 n3=lpby(k+n)
1780 CALL colinear3(n1,n2,n3,x,icol)
1781 IF (icol==0) cycle
1782 n1=n2
1783 n2=n3
1784 END DO
1785 IF (icol==1) dofspc(m)=ndof(m)
1786 END IF
1787 END IF
1788 END IF !(ICOL==1) THEN
1789 ENDDO
1790C------interface 2-----------
1791 DO i=1,nint2
1792 n=iint2(i)
1793 nsn = ipari(5,n)
1794 nmn = ipari(6,n)
1795 ji=ipari(1,n)
1796 k10=ji-1
1797 k11=k10+4*ipari(3,n)
1798C------IRECT(4,NSN)-----
1799 k12=k11+4*ipari(4,n)
1800C------NSV(NSN)--node number---
1801 k13=k12+nsn
1802C------MSR(NMN)-----
1803 k14=k13+nmn
1804C------IRTL(NSN)--main el number---
1805c DO J=1,NSN
1806c NS=INBUF(K12+J)
1807c IF (NDOF(NS)>0) THEN
1808c L=INBUF(K14+J)
1809c ID=K11+4*(L-1)
1810c DO M=1,4
1811c NM=INBUF(ID+M)
1812c DOFSPC(NM)=0
1813c ENDDO
1814c ENDIF
1815c ENDDO
1816 DO j=1,nmn
1817 nm=intbuf_tab(n)%MSR(j)
1818 dofspc(nm)=0
1819 ENDDO
1820 ENDDO
1821C
1822 DO i=1,nrbe3
1823 iad=irbe3(1,i)
1824 ns =irbe3(3,i)
1825 IF (ns==0) cycle
1826 nmn=irbe3(5,i)
1827 DO j=1,nmn
1828 nm=lrbe3(iad+j)
1829 dofspc(nm)=0
1830 ENDDO
1831 ENDDO
1832C
1833 DO i=1,nrbe2
1834 j=irbe2(3,i)
1835 iad=irbe2(1,i)
1836 nsn = irbe2(5,i)
1837 IF (nsn==0) cycle
1838 icol = 1
1839 dofspc(j)=0
1840 DO n =1,nsn
1841 ns = lrbe2(iad+n)
1842 IF (ndof(ns)>3) THEN
1843 icol = 0
1844 cycle
1845 ENDIF
1846 END DO
1847 IF (icol==1) THEN
1848 IF (nsn==1) THEN
1849 dofspc(j)=ndof(j)
1850 ELSEIF (nsn==2) THEN
1851 CALL colinear3(j,lrbe2(iad+1),lrbe2(iad+2),x,icol)
1852 IF (icol==1) dofspc(j)=ndof(j)
1853 ELSE
1854 CALL colinear3(j,lrbe2(iad+1),lrbe2(iad+2),x,icol)
1855 IF (icol==1) THEN
1856 n1=lrbe2(iad+1)
1857 n2=lrbe2(iad+2)
1858 DO n =3,nsn
1859 n3=lrbe2(iad+n)
1860 CALL colinear3(n1,n2,n3,x,icol)
1861 IF (icol==0) cycle
1862 n1=n2
1863 n2=n3
1864 END DO
1865 IF (icol==1) dofspc(j)=ndof(j)
1866 END IF
1867 END IF
1868 END IF !(ICOL==1) THEN
1869 END DO
1870C----6---------------------------------------------------------------7---------8
1871 RETURN
1872 END
1873!||====================================================================
1874!|| colinear3 ../engine/source/implicit/upd_glob_k.F
1875!||--- called by ------------------------------------------------------
1876!|| ini_dofspc ../engine/source/implicit/upd_glob_k.F
1877!||====================================================================
1878 SUBROUTINE colinear3(N1,N2,N3,X,ICOL)
1879C-----------------------------------------------
1880C I m p l i c i t T y p e s
1881C-----------------------------------------------
1882#include "implicit_f.inc"
1883C-----------------------------------------------
1884C C o m m o n B l o c k s
1885C-----------------------------------------------
1886C-----------------------------------------------------------------
1887C D u m m y A r g u m e n t s
1888C-----------------------------------------------
1889 INTEGER N1,N2,N3,ICOL
1890 my_real
1891 . x(3,*)
1892C-----------------------------------------------
1893C L o c a l V a r i a b l e s
1894C-----------------------------------------------
1895 INTEGER I, J ,K
1896 my_real
1897 . x21,y21,z21,x31,y31,z31,pr
1898C
1899 icol=0
1900 x21=x(1,n2)-x(1,n1)
1901 y21=x(2,n2)-x(2,n1)
1902 z21=x(3,n2)-x(3,n1)
1903 x31=x(1,n3)-x(1,n1)
1904 y31=x(2,n3)-x(2,n1)
1905 z31=x(3,n3)-x(3,n1)
1906C -----------------
1907 pr = y21 * z31 - z21 * y31
1908 IF (abs(pr)<=em10) THEN
1909 pr = z21 * x31 - x21 * z31
1910 IF (abs(pr)<=em10) THEN
1911 pr = x21 * y31 - y21 * x31
1912 IF (abs(pr)<=em10) icol=1
1913 END IF
1914 END IF
1915C
1916 RETURN
1917 END
subroutine bc_impa(iadk, jdik, diag_k, lt_k, ndof, iddl, ikc)
Definition bc_imp0.F:1497
subroutine upd_aspc(nddl, ndof, iddl, ikc, itab, iadk, jdik, diag_k, lt_k)
Definition bc_imp0.F:2324
subroutine bc_imp1(icodt, icodr, iskew, skew, ifix, ndof, iadn, iadk, jdik, diag_k, lt_k)
Definition bc_imp0.F:162
subroutine bc_impr1(icodt, icodr, iskew, skew, ndof, iadn, lb)
Definition bc_imp0.F:1124
subroutine bc_imp0(icodt, icodr, iskew, ifix, ndof, iadn)
Definition bc_imp0.F:31
#define my_real
Definition cppsort.cpp:32
if(complex_arithmetic) id
subroutine fv_imp0(iddl, ifix, ndof, iadk, jdik, diag_k, lt_k, ud, nbk, iab, bk, nddl, rd)
Definition fv_imp0.F:37
subroutine dim_fvbcl(ibfv, lj, iskew, icodt, icodr, nddl, iddl, ifix, iadk, jdik, skew, nfvbcl, nbkud)
Definition fv_imp0.F:2300
subroutine fv_rw0(iddl, ifix, ndof, iadk, jdik, diag_k, lt_k, ud, b)
Definition fv_imp0.F:546
subroutine fv_imprl(ibfv, skew, xframe, lj, iddl, ndof, lb)
Definition fv_imp0.F:1555
subroutine fvbc_allo
Definition fv_imp0.F:3203
subroutine fv_impi(iddl, ifix, ndof, iadk, jdik, diag_k, lt_k, ud, b, nddl)
Definition fv_imp0.F:441
subroutine fvbc_impl(ibfv, skew, xframe, lj, iddl, ifix, ndof, iadk, jdik, diag_k, lt_k, ud, rd, lb, nddl, icodt, icodr, icodt1, icodr1, nkud1, ikud, bkud)
Definition fv_imp0.F:2462
subroutine fv_impl(ibfv, skew, xframe, lj, iddl, ifix, ndof, iadk, jdik, diag_k, lt_k, ud, rd, lb)
Definition fv_imp0.F:793
subroutine area(d1, x, x2, y, y2, eint, stif0)
subroutine i2_impi(nint2, iint2, ipari, intbuf_tab, x, ms, in, nss2, iss2, weight, ikc, ndof, nddl, iddl, iadk, jdik, diag_k, lt_k, iaint2, b, itab)
Definition i2_imp0.F:109
subroutine i2_imp0(nint2, iint2, ipari, intbuf_tab, x, ms, in, nmc2, imij2, itab, nsc2, isij2, nss2, iss2, weight, ikc, ndof, nddl, iddl, iadk, jdik, diag_k, lt_k, b)
Definition i2_imp0.F:39
subroutine i2_impr1(ipari, intbuf_tab, x, ndof, iddl, b)
Definition i2_imp1.F:1901
subroutine i2_impr2(ipari, intbuf_tab, a, ar, x, ndof, iddl, b)
Definition i2_imp1.F:2183
subroutine i3cst3(x1, x2, x3, x4, xi, y1, y2, y3, y4, yi, z1, z2, z3, z4, zi, xface, n1, n2, n3, ssc, ttc, x0, y0, z0, xx1, xx2, xx3, xx4, yy1, yy2, yy3, yy4, zz1, zz2, zz3, zz4, xi1, xi2, xi3, xi4, yi1, yi2, yi3, yi4, zi1, zi2, zi3, zi4, xn1, xn2, xn3, xn4, yn1, yn2, yn3, yn4, zn1, zn2, zn3, zn4, area, lft, llt)
Definition i3cst3.F:50
subroutine i3dis3(igimp, nty, dist, x1, x2, x3, x4, xi, y1, y2, y3, y4, yi, z1, z2, z3, z4, zi, xface, n1, n2, n3, ssc, ttc, alp, ans, xp, yp, zp, h1, h2, h3, h4, lft, llt)
Definition i3dis3.F:42
subroutine i3gap3(gap, area, thk, alp, lft, llt)
Definition i3gap3.F:35
subroutine i5cork3(x, irect, msr, nsv, irtl, i_n, i_e, nsvg, jx1, jx2, jx3, jx4, x1, x2, x3, x4, y1, y2, y3, y4, z1, z2, z3, z4, xface, xi, yi, zi, ix1, ix2, ix3, ix4, lft, llt, nft)
Definition i5cork3.F:39
subroutine imp_fhht1(nddl0, nddl, lb, ikc)
Definition imp_dyna.F:2104
subroutine imp_fhht(nddl, lb)
Definition imp_dyna.F:459
subroutine imp_dycrb(am, in, vr, nby, rby0, weight, icodr, iskew, skew)
Definition imp_dyna.F:1211
subroutine imp_dynar(dy_ac, dy_acr, ms, in, fint, mint, v, vr)
Definition imp_dyna.F:393
subroutine imp_setba(f, m, iddl, ndof, b, iflag)
Definition imp_setb.F:135
subroutine imp_setb(f, m, iddl, ndof, b)
Definition imp_setb.F:40
subroutine imp_setbp(f, m, iddl, ndof, ikc, inloc, b)
Definition imp_setb.F:85
subroutine imp_solv(timers, python, icode, iskew, iskwn, ipart, ixtg, ixs, ixq, ixc, ixt, ixp, ixr, ixtg1, itab, itabm1, npc, ibcl, ibfv, sensor_tab, nnlink, lnlink, iparg, igrv, ipari, intbuf_tab, nprw, iconx, npby, lpby, lrivet, nstrf, ljoint, icodt, icodr, isky, adsky, iads_f, ilink, llink, weight, itask, ibvel, lbvel, fbvel, x, d, v, vr, dr, thke, damp, ms, in, pm, skews, geo, eani, bufmat, bufgeo, bufsf, tf, forc, vel, fsav, agrv, fr_wave, parts0, elbuf, rby, rivet, fr_elem, iad_elem, wa, a, ar, stifn, stifr, partsav, fsky, fskyi, iframe, xframe, w16, iactiv, fskym, igeo, ipm, wfext, nodft, nodlt, nint7, num_imp, ns_imp, ne_imp, ind_imp, it, rwbuf, lprw, fr_wall, nbintc, intlist, fopt, rwsav, fsavd, graphe, fac_k, ipiv_k, nkcond, nsensor, monvol, igrsurf, fr_mv, volmon, dirul, nodglob, mumps_par, cddlp, isendto, irecvfrom, newfront, imsch, i2msch, isizxv, ilenxv, islen7, irlen7, islen11, irlen11, islen17, irlen17, irlen7t, islen7t, kinet, num_imp1, temp, dt2prev, waint, lgrav, sh4tree, sh3tree, irlen20, islen20, irlen20t, islen20t, irlen20e, islen20e, irbe3, lrbe3, frbe3, fr_i2m, iad_i2m, fr_rbe3m, iad_rbe3m, frwl6, irbe2, lrbe2, intbuf_tab_c, ikine, diag_sms, icfield, lcfield, cfield, count_remslv, count_remslve, elbuf_tab, elbuf_imp, xdp, weight_md, stack, dimfb, fbsav6, stabsen, tabsensor, drape_sh4n, drape_sh3n, h3d_data, multi_fvm, igrbric, igrsh4n, igrsh3n, igrbeam, forneqs, maxdgap, nddl0, nnzk0, it_t, impbuf_tab, cptreac, fthreac, nodreac, drapeg, interfaces, th_surf, dpl0cld, vel0cld, snpc, stf, glob_therm, wfext_md)
Definition imp_solv.F:173
subroutine imp_chkm(timers, python, icode, iskew, iskwn, ipart, ixtg, ixs, ixq, ixc, ixt, ixp, ixr, ixtg1, itab, itabm1, npc, ibcl, ibfv, sensor_tab, nnlink, lnlink, iparg, igrv, ipari, intbuf_tab, nprw, iconx, npby, lpby, lrivet, nstrf, ljoint, icodt, icodr, isky, adsky, iads_f, ilink, llink, weight, itask, ibvel, lbvel, fbvel, x, d, v, vr, dr, thke, damp, ms, in, pm, skews, geo, eani, bufmat, bufgeo, bufsf, tf, forc, vel, fsav, agrv, fr_wave, parts0, elbuf, rby, rivet, fr_elem, iad_elem, nsensor, wa, a, ar, stifn, stifr, partsav, fsky, fskyi, iframe, xframe, w16, iactiv, fskym, igeo, ipm, wfext, nodft, nodlt, nint7, num_imp, ns_imp, ne_imp, ind_imp, it, rwbuf, lprw, fr_wall, nbintc, intlist, fopt, rwsav, fsavd, dirul, lgrav, irbe3, lrbe3, frbe3, frwl6, irbe2, lrbe2, icfield, lcfield, cfield, elbuf_tab, weight_md, stack, dimfb, fbsav6, stabsen, tabsensor, drape_sh4n, drape_sh3n, h3d_data, nddl0, nnzk0, impbuf_tab, cptreac, fthreac, nodreac, drapeg, th_surf, dpl0cld, vel0cld, snpc, stf, wfext_md, igrsurf)
Definition imp_solv.F:3132
subroutine spmd_sumf_v(v)
Definition imp_spmd.F:1650
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
integer nfvbcl
integer, dimension(:), allocatable icr_1
integer nkud_l
integer, dimension(:), allocatable ikud_1
integer, dimension(:), allocatable ict_1
integer nkud_1
subroutine cp_real(n, x, xc)
Definition produt_v.F:871
subroutine produt_w(nddl, x, y, w, r)
Definition produt_v.F:106
subroutine produt_uhp2(nddl0, nddl, iddl, ndof, ikc, d1, d1r, d2, d2r, norm2, w_imp)
Definition produt_v.F:3398
subroutine produt_v(nddl, x, y, r)
Definition produt_v.F:33
subroutine rbe2_impr1(irbe2, lrbe2, x, skew, ndof, iddl, b, weight)
Definition rbe2_imp0.F:464
subroutine rbe2_imp0(irbe2, lrbe2, x, nsrb2, isb2, ikc, ndof, iddl, iadk, jdik, diag_k, lt_k, b, weight, itab, skew)
Definition rbe2_imp0.F:37
subroutine rbe2_impi(irbe2, lrbe2, x, skew, nsb2, isb2, ikc, ndof, iddl, iadk, jdik, diag_k, lt_k, b, weight, itab)
Definition rbe2_imp0.F:102
subroutine rbe3_imp0(irbe3, lrbe3, frbe3, x, skew, iss3, ikc, ndof, iddl, iadk, jdik, diag_k, lt_k, b, weight, itab)
Definition rbe3_imp0.F:37
subroutine rbe3_impr2(irbe3, lrbe3, frbe3, x, skew, ndof, iddl, b, weight, a, ar)
Definition rbe3_imp0.F:445
subroutine rbe3_impi(irbe3, lrbe3, frbe3, x, skew, nss3, iss3, ikc, ndof, iddl, iadk, jdik, diag_k, lt_k, b, weight, itab)
Definition rbe3_imp0.F:110
subroutine rbe3_impr1(irbe3, lrbe3, frbe3, x, skew, ndof, iddl, b, weight)
Definition rbe3_imp0.F:318
subroutine rby_impr1(x, rby, nod, nby, ndof, iddl, b)
Definition rby_imp0.F:707
subroutine rby_impi(x, rby, lpby, npby, skew, nrbyac, irbyac, nss, iss, iskew, itab, weight, ms, in, nddl, iadk, jdik, diag_k, lt_k, ndof, iddl, ikc, b)
Definition rby_imp0.F:98
subroutine rby_imp0(x, rby, lpby, npby, skew, nrbyac, irbyac, nsc, isij, nmc, imij, nss, iss, iskew, itab, weight, ms, in, nddl, iadk, jdik, diag_k, lt_k, ndof, iddl, ikc, b)
Definition rby_imp0.F:37
subroutine rby_impr2(x, rby, nod, nby, ndof, iddl, b, ac, acr)
Definition rby_imp0.F:820
subroutine recukin(rby, lpby, npby, skew, iskew, itab, weight, ms, in, ibfv, vel, icodt, icodr, nrbyac, irbyac, nint2, iint2, ipari, intbuf_tab, ndof, d, dr, x, xframe, lj, ixr, ixc, ixtg, sh4tree, sh3tree, irbe3, lrbe3, frbe3, irbe2, lrbe2)
Definition recudis.F:103
subroutine rm_imp0(nddl, iadk, jdik, diag_k, lt_k, ndof, iddl, ikc, b, itab)
Definition rm_imp0.F:385
subroutine fv_rwlr0(iddl, b)
Definition srw_imp.F:196
subroutine ud_intg5(lft, llt, nsvg, jx1, jx2, jx3, jx4, ud, h1, h2, h3, h4, xface)
subroutine upd_rhs(icodt, icodr, iskew, ibfv, xframe, rby, x, skew, lpby, npby, nrbyac, irbyac, nint2, iint2, ipari, intbuf_tab, ndof, iddl, ikc, nddl0, b, iupd, inloc, lj, a, ar, ac, acr, nt_rw, irflag, w_ddl, nddl, r02, idyna, v, vr, ms, in, irbe3, lrbe3, frbe3, weight, irbe2, lrbe2)
Definition upd_glob_k.F:664
subroutine ud_int5(ipari, intbuf_tab, x, num_imp, cand_n, cand_e, ud)
subroutine condens_k(nddl, nnz, iadk, jdik, lt_k, ikc)
Definition upd_glob_k.F:222
subroutine upd_int_k(icodt, icodr, iskew, ibfv, npc, tf, vel, xframe, rby, x, skew, lpby, npby, itab, weight, ms, in, nrbyac, irbyac, nss, iss, ipari, intbuf_tab, nint2, iint2, iaint2, nss2, iss2, nddli, nnzi, iadi, jdii, diag_i, lt_i, iddli, nddl, iadk, jdik, ikc, diag_k, lt_k, iddl, ndofi, itok, ud, lb, luj, nt_rw, irbe3, lrbe3, frbe3, nss3, iss3, irbe2, lrbe2, nsb2, isb2)
Definition upd_glob_k.F:465
subroutine condens_ind(nddl, nnz, iadk, jdik, ikc)
Definition upd_glob_k.F:330
subroutine ud_int7(jlt, ns_imp, ne_imp, irect, nsv, nsn, x, ud)
subroutine rer02(rby, lpby, npby, skew, iskew, itab, weight, ms, in, ibfv, vel, icodt, icodr, nrbyac, irbyac, nint2, iint2, ipari, intbuf_tab, ndof, d, dr, x, xframe, lj, ixr, ixc, ixtg, sh4tree, sh3tree, irbe3, lrbe3, frbe3, iadk, jdik, diag_k, lt_k, iddl, ikc, inloc, num_imp, ns_imp, ne_imp, index2, nddl, w_ddl, a, ar, r02, irbe2, lrbe2, x_c)
Definition upd_glob_k.F:940
subroutine colinear3(n1, n2, n3, x, icol)
subroutine ext_rhs(icodt, icodr, iskew, ibfv, xframe, x, skew, ndof, iddl, ikc, nddl0, b, inloc, lj, ac, acr, nt_rw, w_ddl, nddl, r02)
Definition upd_glob_k.F:862
subroutine condens_k0(nddl, nnz, iadk, jdik, nl, lt_k, ikc, ifix, ndf, ndl, nzf, nzl)
Definition upd_glob_k.F:273
subroutine condens_b(nddl, ikc, b)
Definition upd_glob_k.F:400
subroutine inihi11(ns, nm, x, hs, hm)
subroutine inihi7(ns, irect, x, h)
subroutine rer_int_v(ipari, intbuf_tab, num_imp, ns_imp, ne_imp, x, ud)
subroutine upd_glob_k(icodt, icodr, iskew, ibfv, npc, tf, vel, xframe, rby, x, skew, lpby, npby, itab, weight, ms, in, nrbyac, irbyac, nsc, isij, nmc, imij, nss, iss, nint2, iint2, nsc2, isij2, nss2, iss2, ipari, intbuf_tab, nddl, nnz, iadk, jdik, diag_k, lt_k, ndof, iddl, ikc, ud, b, nkud, ikud, bkud, nmc2, imij2, nt_rw, rd, lj, irbe3, lrbe3, frbe3, iss3, irbe2, lrbe2, isb2, nsrb2)
Definition upd_glob_k.F:66
subroutine ud_int11(jlt, ns_imp, ne_imp, irects, irectm, nsn, x, ud)
subroutine ini_dofspc(npby, lpby, nrbyac, irbyac, nint2, iint2, ipari, intbuf_tab, ndof, irbe3, lrbe3, irbe2, lrbe2, x, dofspc)