OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
hm_read_rbe3.F File Reference
#include "implicit_f.inc"
#include "scr17_c.inc"
#include "com01_c.inc"
#include "com04_c.inc"
#include "units_c.inc"
#include "param_c.inc"
#include "tabsiz_c.inc"
#include "r2r_c.inc"
#include "sphcom.inc"
#include "scr03_c.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine hm_read_rbe3 (irbe3, lrbe3, frbe3, itab, itabm1, igrnod, iskn, lxintd, ikine, iddlevel, nom_opt, itagnd, grnod_uid, unitab, lsubmodel)
subroutine inirbe3 (irbe3, lrbe3, frbe3, skew, x, ms, in, nom_opt, stifn, stifr, icode, itab)
subroutine rbe3chk (inrbe3, ilrbe3, ns, xyz, frbe3, skew, ng, irot, imodif, wmin, ipen, ierr)
subroutine rbe3uf (inrbe3, ilrbe3, el, tw, xyz, refpt, fufxlc, fufylc, fufzlc, fufx, fufy, fufz, mufx, mufy, mufz, tfufx, tfufy, tfufz, tmufx, tmufy, tmufz, denfx, denfy, denfz, ng)
subroutine rbe3um (inrbe3, ilrbe3, el, tw, rw, xyz, refpt, cgmx, cgmy, cgmz, fumxlc, fumylc, fumzlc, mxlc, mylc, mzlc, fumx, fumy, fumz, mx, my, mz, mumx, mumy, mumz, tfumx, tfumy, tfumz, tmumx, tmumy, tmumz, averef, denmx, denmy, denmz, ng, irot)
subroutine invert (matrix, inverse, n, errorflag)
subroutine wrrinf (title, r, n)
subroutine zero1 (r, n)
subroutine prerbe3fr (irbe3, n, jt, jr)
subroutine hireorbe3 (irbe3, lrbe3, frbe3, nom_opt)

Function/Subroutine Documentation

◆ hireorbe3()

subroutine hireorbe3 ( integer, dimension(nrbe3l,*) irbe3,
integer, dimension(*) lrbe3,
frbe3,
integer, dimension(lnopt1,*) nom_opt )

Definition at line 1593 of file hm_read_rbe3.F.

1594C-----------------------------------------------
1595C M o d u l e s
1596C-----------------------------------------------
1597 USE my_alloc_mod
1598 USE message_mod
1599 USE names_and_titles_mod , ONLY : nchartitle
1600C-----------------------------------------------
1601C I m p l i c i t T y p e s
1602C-----------------------------------------------
1603#include "implicit_f.inc"
1604C-----------------------------------------------
1605C C o m m o n B l o c k s
1606C-----------------------------------------------
1607#include "com04_c.inc"
1608#include "param_c.inc"
1609#include "scr17_c.inc"
1610#include "tabsiz_c.inc"
1611C-----------------------------------------------
1612C D u m m y A r g u m e n t s
1613C-----------------------------------------------
1614 INTEGER IRBE3(NRBE3L,*), LRBE3(*),NOM_OPT(LNOPT1,*)
1615 my_real
1616 . frbe3(*)
1617C-----------------------------------------------
1618C L o c a l V a r i a b l e s
1619C-----------------------------------------------
1620 INTEGER I, N, J,K, NS,NM, NI, NMT,M,NFT,IAD,II,M2,III,IPEN
1621C
1622 INTEGER ID, INDEX(NRBE3),NZ,IAD1,IADS,
1623 . LCOPY(SLRBE3),ICOPY(NRBE3L,NRBE3)
1624 INTEGER,DIMENSION(:),ALLOCATABLE :: ITAG
1625 CHARACTER(LEN=NCHARTITLE)::TITR
1626 my_real fcopy(slrbe3*3)
1627C========================================================================|
1628C-----only one level of hierarchy is allowed
1629 CALL my_alloc(itag,numnod)
1630 itag(1:numnod)=0
1631 DO i=1,nrbe3
1632 ns = irbe3(3,i)
1633 IF (itag(ns)==0) itag(ns)=i
1634 index(i) = i
1635 ENDDO
1636 nz = 0
1637! case Iform=2 : error out if >1 level
1638 DO i=1,nrbe3
1639 iad = irbe3(1,i)
1640 nm = irbe3(5,i)
1641 ipen= irbe3(9,i)
1642 IF (ipen>=0) cycle
1643 DO j =1,nm
1644 m = lrbe3(iad+j)
1645 IF (itag(m)>0) THEN
1646 ii = itag(m)
1647 IF (irbe3(9,ii)<0) nz = nz +1
1648 END IF
1649 ENDDO
1650 ENDDO
1651 IF (nz >0 ) THEN
1652C-----error out if >1 level
1653 DO i=1,nrbe3
1654 iad = irbe3(1,i)
1655 nm = irbe3(5,i)
1656 ipen= irbe3(9,i)
1657 IF (ipen>=0) cycle
1658 DO j =1,nm
1659 m = lrbe3(iad+j)
1660 IF (itag(m)>0) THEN
1661 ii = itag(m)
1662 IF (irbe3(9,ii)<0) THEN
1663 nmt = irbe3(5,ii)
1664 DO k =1,nmt
1665 nft = irbe3(1,ii)
1666 m2 = lrbe3(nft+k)
1667 IF (itag(m2)>0) THEN
1668 iii = itag(m2)
1669 IF (irbe3(9,iii)<0) THEN
1670 id=nom_opt(1,iii)
1671 CALL fretitl2(titr,nom_opt(lnopt1-ltitr+1,iii),ltitr)
1672 CALL ancmsg(msgid=1887,
1673 . msgtype=msgerror,
1674 . anmode=aninfo,
1675 . i1=id,
1676 . c1=titr,
1677 . i2=nom_opt(1,i),
1678 . i3=nom_opt(1,ii))
1679 END IF !(IRBE3(9,III)<0) THEN
1680 END IF
1681 ENDDO
1682 END IF !(IRBE3(9,II)<0) THEN
1683 ENDIF
1684 ENDDO
1685 ENDDO
1686 END IF !(NZ >0 ) THEN
1687! case Iform=1,3
1688 nz = 0
1689 DO i=1,nrbe3
1690 iad = irbe3(1,i)
1691 ipen= irbe3(9,i)
1692 IF (ipen>0) cycle
1693 nm = irbe3(5,i)
1694 DO j =1,nm
1695 m = lrbe3(iad+j)
1696 IF (itag(m)>0) THEN
1697 ipen= irbe3(9,itag(m))
1698 IF (ipen==0) THEN
1699 nz = nz + 1
1700 cycle
1701 END IF
1702 ENDIF
1703 ENDDO
1704 ENDDO
1705 IF (nz >0 ) THEN
1706!-----swtch to penalty if >1 level
1707 DO i=nrbe3,1,-1
1708 iad = irbe3(1,i)
1709 ipen= irbe3(9,i)
1710 IF (ipen>0) cycle
1711 nm = irbe3(5,i)
1712 ii = 0
1713 DO j =1,nm
1714 m = lrbe3(iad+j)
1715 IF (itag(m)>0) THEN
1716 ipen= irbe3(9,itag(m))
1717 IF (ipen/=0) cycle
1718 ii = itag(m)
1719 nmt = irbe3(5,ii)
1720 nft = irbe3(1,ii)
1721 iii = 0
1722 DO k =1,nmt
1723 m2 = lrbe3(nft+k)
1724 IF (itag(m2)>0) THEN
1725 ipen= irbe3(9,itag(m2))
1726 IF (ipen==0) THEN
1727 iii = itag(m2)
1728 cycle
1729 END IF
1730 END IF
1731 ENDDO
1732 IF (iii>0) cycle
1733 ENDIF
1734 ENDDO
1735 IF (ii >0 .AND. iii>0) THEN
1736 id=nom_opt(1,i)
1737 CALL fretitl2(titr,nom_opt(lnopt1-ltitr+1,i),ltitr)
1738 irbe3(9,i) = 1
1739 CALL ancmsg(msgid=3099,
1740 . msgtype=msgwarning,
1741 . anmode=aninfo,
1742 . i1=id,
1743 . c1=titr,
1744 . i2=nom_opt(1,ii),
1745 . i3=nom_opt(1,iii))
1746 END IF
1747 ENDDO
1748 END IF
1749C----- re-ordering
1750 nz = 0
1751 DO i=1,nrbe3
1752 iad = irbe3(1,i)
1753 ipen= irbe3(9,i)
1754 IF (ipen>0) cycle
1755 nm = irbe3(5,i)
1756 DO j =1,nm
1757 m = lrbe3(iad+j)
1758 IF (itag(m)>0) THEN
1759 ipen= irbe3(9,itag(m))
1760 IF (ipen<=0) THEN
1761 nz = nz +1
1762!----- exchange INDEX(I) , NZ
1763 IF (index(i) > nz) THEN
1764 ni = index(i)
1765 index(i) = nz
1766 index(nz) = ni
1767 END IF
1768 cycle
1769 END IF !(IPEN==0) THEN
1770 ENDIF
1771 ENDDO
1772 ENDDO
1773 IF (nz >0 ) THEN
1774 iads = slrbe3/2
1775 lcopy(1:slrbe3) = lrbe3(1:slrbe3)
1776 icopy(1:nrbe3l,1:nrbe3) = irbe3(1:nrbe3l,1:nrbe3)
1777 fcopy(1:slrbe3*3) = frbe3(1:slrbe3*3)
1778 iad1 = 0
1779 DO n=1,nrbe3
1780 i = index(n)
1781 iad = icopy(1,i)
1782 ns = icopy(3,i)
1783 nm = icopy(5,i)
1784 irbe3(1,n) = iad1
1785 DO j =2,nrbe3l
1786 irbe3(j,n) = icopy(j,i)
1787 ENDDO
1788 nom_opt(1,n)=irbe3(2,n)
1789 DO j =1,nm
1790 lrbe3(iad1+j)=lcopy(iad+j)
1791 lrbe3(iads+iad1+j)=lcopy(iads+iad+j)
1792 ENDDO
1793 DO j =1,6*nm
1794 frbe3(6*iad1+j)=fcopy(6*iad+j)
1795 ENDDO
1796 iad1 =iad1+nm
1797 ENDDO
1798 END IF
1799C
1800 DEALLOCATE(itag)
1801 RETURN
#define my_real
Definition cppsort.cpp:32
subroutine icopy(n, sx, incx, sy, incy)
ICOPY
Definition icopy.f:75
initmumps id
integer, parameter nchartitle
subroutine ancmsg(msgid, msgtype, anmode, i1, i2, i3, i4, i5, i6, i7, i8, i9, i10, i11, i12, i13, i14, i15, i16, i17, i18, i19, i20, r1, r2, r3, r4, r5, r6, r7, r8, r9, c1, c2, c3, c4, c5, c6, c7, c8, c9, prmode)
Definition message.F:895
subroutine fretitl2(titr, iasc, l)
Definition freform.F:799

◆ hm_read_rbe3()

subroutine hm_read_rbe3 ( integer, dimension(nrbe3l,*) irbe3,
integer, dimension(*) lrbe3,
frbe3,
integer, dimension(*) itab,
integer, dimension(*) itabm1,
type (group_), dimension(ngrnod) igrnod,
integer, dimension(liskn,*) iskn,
integer lxintd,
integer, dimension(*) ikine,
integer iddlevel,
integer, dimension(lnopt1,*) nom_opt,
integer, dimension(*) itagnd,
integer grnod_uid,
type (unit_type_), intent(in) unitab,
type(submodel_data), dimension(*), intent(in) lsubmodel )

Definition at line 45 of file hm_read_rbe3.F.

48C-------------------------------------
49C READING STRUCTURE RIGIDES
50C-----------------------------------------------
51C M o d u l e s
52C-----------------------------------------------
53 USE my_alloc_mod
54 USE unitab_mod
55 USE r2r_mod
56 USE message_mod
57 USE groupdef_mod
58 USE submodel_mod
61C-----------------------------------------------
62C I m p l i c i t T y p e s
63C-----------------------------------------------
64#include "implicit_f.inc"
65C-----------------------------------------------
66C C o m m o n B l o c k s
67C-----------------------------------------------
68#include "scr17_c.inc"
69#include "com01_c.inc"
70#include "com04_c.inc"
71#include "units_c.inc"
72#include "param_c.inc"
73#include "tabsiz_c.inc"
74#include "r2r_c.inc"
75#include "sphcom.inc"
76#include "scr03_c.inc"
77C-----------------------------------------------
78C D u m m y A r g u m e n t s
79C-----------------------------------------------
80 INTEGER IRBE3(NRBE3L,*), LRBE3(*), ITAB(*),ITABM1(*),
81 . ISKN(LISKN,*),LXINTD,
82 . IDDLEVEL,IKINE(*),ITAGND(*)
83 TYPE (UNIT_TYPE_),INTENT(IN) ::UNITAB
84 my_real frbe3(6,*)
85 INTEGER NOM_OPT(LNOPT1,*)
86C-----------------------------------------------
87 TYPE (GROUP_) , DIMENSION(NGRNOD) :: IGRNOD
88 INTEGER :: GRNOD_UID
89 TYPE(SUBMODEL_DATA),INTENT(IN)::LSUBMODEL(*)
90C-----------------------------------------------
91C L o c a l V a r i a b l e s
92C-----------------------------------------------
93 INTEGER I, N, K, NSL, NUSER, NM, NI, NI_OK,
94 . ISK, INGU, IGM, J, IAD,NS,NN,J6(6),JJ,II,
95 . IC,IC1,IC2,IROT,ISKS,IADS,IMODIF,
96 . IDIR,NRB,ID,UID,SUB_INDEX,IFORM
97 my_real w
98 INTEGER, DIMENSION(:),ALLOCATABLE :: IKINE1
99 INTEGER, DIMENSION(:),ALLOCATABLE :: ISKEW0
100 my_real, DIMENSION(:,:),ALLOCATABLE :: wi
101 CHARACTER(LEN=NCHARTITLE) :: TITR
102 CHARACTER(LEN=NCHARFIELD) :: STRING
103 CHARACTER :: MESS*40
104 LOGICAL IS_AVAILABLE
105C-----------------------------------------------
106C E x t e r n a l F u n c t i o n s
107C-----------------------------------------------
108 INTEGER USR2SYS
109C
110 DATA mess/'INTERPOLATION CONSTRAINT BODY '/
111C-----------------------------------------------
112C IRBE3(1,I) : IAD0 for LRBE3 and FRBE3
113C IRBE3(2,I) : TYPE usr' id temporaire (print)
114C IRBE3(3,I) : SECONDARY NODES
115C IRBE3(4,I) : REF_DOF
116C IRBE3(5,I) : NUMBER MAIN NODE
117C IRBE3(6,I) : IROT =0 no rotational dependent
118C IRBE3(7,I) : SENSOR NUMBER - not used yet
119C IRBE3(8,I) : Imodif ! 4 : switch to penalty formulation (for testing)
120C IRBE3(9,I) : if penalty formulation : 0=kinematic, 1=penalty
121C IRBE3(10,I) : Internal number of IRBE3 (necessary for Modif/SPMD)
122C=======================================================================
123! complete case of switching to penalty formulation; free nodes w/ spc?
124! 1) hierachy level>1
125! 2) incompatibility kin
126! 3) avoid too much added mass : 2 conditions 1)Is>0 2)rR>?
127! rbe3_mod is used from engine, ini is done in resol_init,internal arrays of pen are added in restart
128 CALL my_alloc(ikine1,3*numnod)
129 CALL my_alloc(iskew0,slrbe3/2)
130 CALL my_alloc(wi,6,numnod)
131c
132 WRITE(iout,1000)
133CC
134 is_available = .false.
135 CALL hm_option_start('/RBE3')
136C
137 DO i=1,3*numnod
138 ikine1(i) = 0
139 ENDDO
140C
141 iads = slrbe3/2
142 iad = 0
143 nrb = 0
144
145C---------otherwise quadratic cycling---------
146 DO j=1,6
147 DO n=1,numnod
148 wi(j,n)=zero
149 ENDDO
150 ENDDO
151
152 DO i=1,nrbe3
153 nrb=nrb+1
154C----------Multidomaines --> on ignore les rbe3 non tages---------
155 IF(nsubdom>0)THEN
156 IF(tagrb3(nrb)==0)CALL hm_sz_r2r(tagrb3,nrb,lsubmodel)
157 END IF
158 iform = 0
159C--------------------------------------------------
160C EXTRACT DATAS OF /RBE3/... LINE
161C--------------------------------------------------
162 CALL hm_option_read_key(lsubmodel,
163 . option_id = id,
164 . unit_id = uid,
165 . option_titr = titr,
166 . submodel_index = sub_index)
167 nom_opt(1,i)=id
168 CALL fretitl(titr,nom_opt(lnopt1-ltitr+1,i),ltitr)
169C
170 CALL hm_get_intv('dependentnode',nsl,is_available,lsubmodel)
171 CALL hm_get_intv('LTX',j6(1),is_available,lsubmodel)
172 CALL hm_get_intv('LTY',j6(2),is_available,lsubmodel)
173 CALL hm_get_intv('LTZ',j6(3),is_available,lsubmodel)
174 CALL hm_get_intv('LRX',j6(4),is_available,lsubmodel)
175 CALL hm_get_intv('LRY',j6(5),is_available,lsubmodel)
176 CALL hm_get_intv('LRZ',j6(6),is_available,lsubmodel)
177 CALL hm_get_intv('nset',nn,is_available,lsubmodel)
178 CALL hm_get_intv('I_Modif',imodif,is_available,lsubmodel)
179 CALL hm_get_intv('Iform',iform,is_available,lsubmodel)
180 irbe3(2,i) = id
181 irbe3(10,i)=i
182 nuser = id
183C
184 IF (imodif==0) imodif =1
185! IF (IMODIF==4) IRBE3(9,I) =1
186 irbe3(8,i) = imodif
187 SELECT CASE (iform)
188 CASE(0,1)
189 irbe3(9,i) =0
190 IF (iform==0) iform =1
191 CASE(2)
192 irbe3(9,i) = -1
193 CASE(3)
194 irbe3(9,i) = 1
195 END SELECT
196 ns = usr2sys(nsl,itabm1,mess,nuser)
197 ic1=j6(1)*4 +j6(2)*2 +j6(3)
198 ic2=j6(4)*4 +j6(5)*2 +j6(6)
199 ic =ic1*512+ic2*64
200 IF (ic==0) ic =7*512+7*64
201 IF (ns10e > 0) THEN
202 IF(itagnd(ns)/=0) THEN
203C------- error out
204 CALL ancmsg(msgid=1208,
205 . msgtype=msgerror,
206 . anmode=aninfo_blind_1,
207 . i1=itab(ns),
208 . c1='RBE3 ',
209 . i2=nuser,
210 . c2='RBE3 ')
211 ENDIF
212 END IF
213 irbe3(3,i) = ns
214 irbe3(4,i) = ic
215 irbe3(1,i) = iad
216 irot = 0
217 DO j=1,nn
218 CALL hm_get_float_array_index('independentnodesetcoeffs',w,j,is_available,lsubmodel,unitab)
219 CALL hm_get_int_array_index('tx',j6(1),j,is_available,lsubmodel)
220 CALL hm_get_int_array_index('ty',j6(2),j,is_available,lsubmodel)
221 CALL hm_get_int_array_index('tz',j6(3),j,is_available,lsubmodel)
222 CALL hm_get_int_array_index('rx',j6(4),j,is_available,lsubmodel)
223 CALL hm_get_int_array_index('ry',j6(5),j,is_available,lsubmodel)
224 CALL hm_get_int_array_index('rz',j6(6),j,is_available,lsubmodel)
225 CALL hm_get_int_array_index('SKEW_ARRAY',isk,j,is_available,lsubmodel)
226 CALL hm_get_int_array_index('independentnodesets',ingu,j,is_available,lsubmodel)
227 IF (w==zero.OR.imodif==3) w=one
228 nm = 0
229 IF(ingu > 0) THEN
230 CALL c_hash_find(grnod_uid,ingu,igm)
231 IF (igm == 0)THEN
232 CALL ancmsg(msgid=53,
233 . msgtype=msgerror,
234 . anmode=aninfo,
235 . c1= mess,
236 . i1=ingu)
237 ELSE
238 nm = igrnod(igm)%NENTITY
239 lrbe3(iad+1:iad+nm) = igrnod(igm)%ENTITY(1:nm)
240 ENDIF
241 END IF !(INGU > 0)
242
243 isks = 0
244 IF ((j6(1)+j6(2)+j6(3)+j6(4)+j6(5)+j6(6))==0) THEN
245 j6(1)=1
246 j6(2)=1
247 j6(3)=1
248 ENDIF
249 IF (isk/=0) THEN
250 DO jj=0,numskw+min(1,nspcond)*numsph+nsubmod
251 IF(isk==iskn(4,jj+1)) THEN
252 isks=jj+1
253 GO TO 10
254 ENDIF
255 ENDDO
256 CALL ancmsg(msgid=184,
257 . msgtype=msgerror,
258 . anmode=aninfo,
259 . c1='RBE3',
260 . i1=nuser,
261 . c2='RBE3',
262 . c3=titr,
263 . i2=isk)
264 10 CONTINUE
265 ENDIF
266 DO k=1,nm
267 ni=lrbe3(iad+k)
268 lrbe3(iads+iad+k)=isks
269 iskew0(iad+k)=isk
270 DO jj=1,6
271 ii = j6(jj)
272 IF (ii>zero) THEN
273 IF (jj>3) irot=1
274 IF (wi(jj,ni)==zero) THEN
275 wi(jj,ni) = w
276 ELSE
277 CALL ancmsg(msgid=705,
278 . msgtype=msgerror,
279 . anmode=aninfo,
280 . i1=nuser,
281 . c1=titr,
282 . i2=itab(ni))
283 ENDIF
284 ENDIF
285 IF (ns10e > 0) THEN
286 IF(itagnd(ni)/=0) THEN
287C------- error out
288 CALL ancmsg(msgid=1211,
289 . msgtype=msgerror,
290 . anmode=aninfo,
291 . i1=itab(ni),
292 . c1='RBE3',
293 . i2=nuser,
294 . c2='RBE3')
295 ENDIF
296 END IF
297 frbe3(jj,iad+k) = wi(jj,ni)
298 ENDDO
299 END DO ! K=1,NM
300 iad = iad+nm
301 ENDDO ! DO J=1,NN
302 irbe3(5,i) = iad-irbe3(1,i)
303 irbe3(6,i) = irot
304
305! optimisation: only concerned nodes are set to zero
306! avoid to set to 0 all nodes (quadratic loop on NRBE3 and NUMNOD)
307 DO ni_ok = irbe3(1,i)+1,irbe3(1,i)+irbe3(5,i)
308 DO jj = 1,6
309 wi(jj,lrbe3(ni_ok)) = zero
310 ENDDO
311 ENDDO
312
313 END DO !I=1,NRBE3
314C
315 IF (ipri<5) WRITE(iout,1103)
316 DO i=1,nrbe3
317 iad = irbe3(1,i)
318 ns = irbe3(3,i)
319 nm = irbe3(5,i)
320 nuser = irbe3(2,i)
321 imodif = irbe3(8,i)
322 IF (imodif/=2) irbe3(8,i)=4
323 CALL prerbe3fr(irbe3 ,i ,j6 ,j6(4) )
324 IF (ipri>=5) THEN
325 WRITE(iout,1100) nuser,itab(ns),j6,nm,imodif,iform
326 WRITE(iout,1101)
327 DO j = 1, nm
328 WRITE(iout,1102) itab(lrbe3(iad+j)),iskew0(iad+j),
329 . (frbe3(jj,iad+j),jj=1,6)
330 ENDDO
331 WRITE(iout,*)
332 ELSE
333c WRITE(IOUT,1103)
334 WRITE(iout,1104) nuser,itab(ns),j6,nm,imodif,iform
335 END IF
336 lxintd = lxintd + nm/4 + 1
337 IF (iddlevel == 0) THEN
338 DO idir=1,6
339 CALL kinset(4096,itab(ns),ikine(ns),idir,0,
340 . ikine1(ns))
341 ENDDO
342 ENDIF
343 ENDDO
344 IF (nspmd==1) lxintd = 0
345C
346 DEALLOCATE(ikine1)
347 DEALLOCATE(iskew0)
348 DEALLOCATE(wi)
349 RETURN
350C
351 1000 FORMAT(//
352 .' INTERPOLATION CONSTRAINT BODY (RBE3) '/
353 . ' ---------------------- '/)
354 1100 FORMAT( 5x,'NUMBER. . . . . . . . . . . . . ',i10
355 . /5x,'DEPENDENT NODE . . . . . . . . .',i10
356 . /5x,'REFERENCE DOF(Trarot). . . . . . . ',3i1,1x,3i1
357 . /5x,'NUMBER OF INDEPENDENT NODES. . .',i10
358 . /5x,'FLAG OF WEIGHTING MODIFICATION .',i10
359 . /5x,'FLAG OF RBE3 FORMULATION. . . . ',i10)
360 1101 FORMAT(//
361 .' WEIGHTING FACTORS OF INDEPENDENT NODES '/
362 .' ------------------- '/
363 .' NODE SKEW DIR_TRA_1 DIR_TRA_2',
364 .' DIR_TRA_3 DIR_ROT_1 DIR_ROT_2',
365 .' DIR_ROT_3'/)
366 1102 FORMAT(3x,2i10,3x,6g20.13)
367 1103 FORMAT(' RBE3_ID DEPENDENT_NODE REF_DOF #IND. IMODIF IFORM'/)
368 1104 FORMAT(3x,2i10,2x,3i1,1x,3i1,3i10)
369C
void c_hash_find(int *map, int *key, int *val)
subroutine hm_get_float_array_index(name, rval, index, is_available, lsubmodel, unitab)
subroutine hm_get_int_array_index(name, ival, index, is_available, lsubmodel)
subroutine hm_get_intv(name, ival, is_available, lsubmodel)
subroutine hm_option_start(entity_type)
subroutine prerbe3fr(irbe3, n, jt, jr)
subroutine kinset(ik, node, ikine, idir, isk, ikine1)
Definition kinset.F:57
#define min(a, b)
Definition macros.h:20
integer, parameter ncharkey
integer, parameter ncharfield
integer, dimension(:), allocatable tagrb3
Definition r2r_mod.F:138
integer nsubmod
subroutine hm_sz_r2r(tag, val, lsubmodel)
subroutine fretitl(titr, iasc, l)
Definition freform.F:615
integer function usr2sys(iu, itabm1, mess, id)
Definition sysfus.F:146

◆ inirbe3()

subroutine inirbe3 ( integer, dimension(nrbe3l,*) irbe3,
integer, dimension(*) lrbe3,
frbe3,
skew,
x,
ms,
in,
integer, dimension(lnopt1,*) nom_opt,
intent(in) stifn,
intent(in) stifr,
integer, dimension(numnod), intent(in) icode,
integer, dimension(numnod), intent(in) itab )

Definition at line 385 of file hm_read_rbe3.F.

388C-----------------------------------------------
389C M o d u l e s
390C-----------------------------------------------
391 USE message_mod
393C-----------------------------------------------
394C I m p l i c i t T y p e s
395C-----------------------------------------------
396#include "implicit_f.inc"
397C-----------------------------------------------
398C C o m m o n B l o c k s
399C-----------------------------------------------
400#include "scr17_c.inc"
401#include "com01_c.inc"
402#include "com04_c.inc"
403#include "param_c.inc"
404#include "tabsiz_c.inc"
405C-----------------------------------------------
406C D u m m y A r g u m e n t s
407C-----------------------------------------------
408 INTEGER IRBE3(NRBE3L,*), LRBE3(*)
409 my_real
410 . skew(lskew,*), frbe3(*),x(3,*),ms(*),in(*)
411 my_real , DIMENSION(NUMNOD), INTENT(IN ) ::stifn
412 my_real , DIMENSION(NUMNOD), INTENT(IN ) ::stifr
413 INTEGER , DIMENSION(NUMNOD), INTENT(IN ) ::ICODE
414 INTEGER , DIMENSION(NUMNOD), INTENT(IN ) ::ITAB
415 INTEGER NOM_OPT(LNOPT1,*)
416C-----------------------------------------------
417C F u n c t i o n
418C-----------------------------------------------
419 INTEGER NLOCAL
420 EXTERNAL nlocal
421C-----------------------------------------------
422C L o c a l V a r i a b l e s
423C-----------------------------------------------
424 INTEGER I, NM, NMT,M,IROT,IMO,ICR,
425 . J, P,IAD,NS,IADS,IERR1,IC,ICT
426 INTEGER ID,IPEN
427 CHARACTER(LEN=NCHARTITLE)::TITR
428C
429 my_real wmin,in_t
430C-----------------------------------------------
431C FRBE3(1-3,I) : tw(i)
432C FRBE3(4-6,I) : rw(i)
433C FRBE3(3*SLRBE3+I) : MS
434C FRBE3(3*SLRBE3+SLRBE3/2+I) : IN
435C========================================================================|
436C----- re-ordering by hierrrchy (only 1 level)
437 CALL hireorbe3(irbe3 ,lrbe3 ,frbe3 ,nom_opt)
438 nmt = slrbe3/2
439 iads = 3*slrbe3
440 DO i=1,nrbe3
441 iad = irbe3(1,i)
442 ns = irbe3(3,i)
443 nm = irbe3(5,i)
444 irot =irbe3(6,i)
445 imo = irbe3(8,i)
446 ipen = irbe3(9,i)
447
448 id=nom_opt(1,i)
449 CALL fretitl2(titr,nom_opt(lnopt1-ltitr+1,i),ltitr)
450! Switch to penalty if /BCS of NS
451 IF (icode(ns)>0 .AND. ipen==0) THEN
452 irbe3(9,i) = 1 ! add message
453 CALL ancmsg(msgid=3115,
454 . msgtype=msgwarning,
455 . anmode=aninfo_blind_1,
456 . i1=id,
457 . c1=titr,
458 . i2=itab(ns))
459 END IF
460C------ check if solid node for all independent
461! will be switch to pen if arm>0
462 IF (in(ns)>em10.AND.min(stifn(ns),stifr(ns))>em10.AND.ipen==0) ipen = 2
463 IF (irot>0) THEN
464 in_t = zero
465 DO j=1,nm
466 m= lrbe3(iad+j)
467 in_t = in_t + in(m)
468 ENDDO
469 IF (in_t<em20) THEN
470 CALL ancmsg(msgid=3009,
471 . msgtype=msgwarning,
472 . anmode=aninfo_blind_1,
473 . i1=id,
474 . c1=titr)
475 irot = 0
476 irbe3(6,i) = 0
477 END IF
478 END IF !(IROT>0) THEN
479 CALL rbe3chk(lrbe3(iad+1),lrbe3(nmt+iad+1) ,ns ,x ,
480 . frbe3(6*iad+1),skew ,nm ,irot ,imo, wmin,
481 . ipen, ierr1 )
482 IF (ipen==-2) THEN
483 irbe3(9,i) = 0 ! small arm keep kinematic
484 ELSEIF (ipen==2) THEN
485 irbe3(9,i) = 1 ! add message
486 CALL ancmsg(msgid=3097,
487 . msgtype=msgwarning,
488 . anmode=aninfo_blind_1,
489 . i1=id,
490 . c1=titr)
491 END IF
492 IF (ierr1==400) THEN
493 id= irbe3(2,i)
494 CALL ancmsg(msgid=3098,
495 . msgtype=msgerror,
496 . anmode=aninfo,
497 . i1=id,
498 . c1=titr)
499 ELSEIF (ierr1>0) THEN
500 id= irbe3(2,i)
501 CALL ancmsg(msgid=706,
502 . msgtype=msgerror,
503 . anmode=aninfo,
504 . i1=id,
505 . c1=titr)
506 ELSEIF (ierr1<0) THEN
507 IF (imo==2) THEN
508 id= irbe3(2,i)
509 CALL ancmsg(msgid=749,
510 . msgtype=msgwarning,
511 . anmode=aninfo_blind_1,
512 . i1=id,
513 . c1=titr)
514 ELSE
515 id= irbe3(2,i)
516 CALL ancmsg(msgid=757,
517 . msgtype=msgwarning,
518 . anmode=aninfo_blind_2,
519 . i1=id,
520 . c1=titr,
521 . r1=wmin)
522 END IF
523 END IF
524 DO j=1,nm
525 m= lrbe3(iad+j)
526 frbe3(iads+j) =ms(m)
527 IF (iroddl/=0) frbe3(iads+nmt+j) =in(m)
528 ENDDO
529C------IF NS not in any proc : done already in domdec2
530 DO p = 1, nspmd
531 IF (nlocal(ns,p)==0) THEN
532 GO TO 200
533 ENDIF
534 DO j=1,nm
535 m= lrbe3(iad+j)
536 IF (nspmd>1) CALL ifrontplus(m,p)
537 ENDDO
538C optimisation possible
539 200 CONTINUE
540 END DO
541 iads = iads + nm
542 ENDDO
543C
544 RETURN
integer function nlocal(n, p)
Definition ddtools.F:350
subroutine ifrontplus(n, p)
Definition frontplus.F:101
subroutine rbe3chk(inrbe3, ilrbe3, ns, xyz, frbe3, skew, ng, irot, imodif, wmin, ipen, ierr)
subroutine hireorbe3(irbe3, lrbe3, frbe3, nom_opt)

◆ invert()

subroutine invert ( dimension(n,n) matrix,
dimension(n,n) inverse,
integer, intent(in) n,
integer, intent(out) errorflag )

Definition at line 1367 of file hm_read_rbe3.F.

1368C-----------------------------------------------
1369C I m p l i c i t T y p e s
1370C-----------------------------------------------
1371#include "implicit_f.inc"
1372c !DECLARATIONS
1373 INTEGER, INTENT(IN) :: N
1374 INTEGER, INTENT(OUT) :: ERRORFLAG !RETURN ERROR STATUS. -1 FOR ERROR, 0 FOR NORMAL
1375 my_real
1376 * , INTENT(IN), DIMENSION(N,N) :: matrix !INPUT MATRIX
1377 my_real
1378 * , INTENT(OUT), DIMENSION(N,N) :: inverse !INVERTED MATRIX
1379
1380 LOGICAL :: FLAG = .true.
1381 INTEGER :: I, J, K
1382 my_real
1383 * :: m
1384 my_real
1385 * , DIMENSION(N,2*N) :: augmatrix !AUGMENTED MATRIX
1386
1387c !AUGMENT INPUT MATRIX WITH AN IDENTITY MATRIX
1388 DO i = 1, n
1389 DO j = 1, 2*n
1390 IF (j <= n ) THEN
1391 augmatrix(i,j) = matrix(i,j)
1392 ELSE IF ((i+n) == j) THEN
1393 augmatrix(i,j) = one
1394 ELSE
1395 augmatrix(i,j) = zero
1396 ENDIF
1397 END DO
1398 END DO
1399
1400c !REDUCE AUGMENTED MATRIX TO UPPER TRIANGULAR FORM
1401 DO k =1, n-1
1402 IF (augmatrix(k,k) == 0) THEN
1403 flag = .false.
1404 DO i = k+1, n
1405 IF (augmatrix(i,k) /= 0) THEN
1406 DO j = 1,2*n
1407 augmatrix(k,j) = augmatrix(k,j)+augmatrix(i,j)
1408 END DO
1409 flag = .true.
1410 EXIT
1411 ENDIF
1412 IF (flag .EQV. .false.) THEN
1413 print*, "MATRIX IS NON - INVERTIBLE"
1414 inverse = 0
1415 errorflag = -1
1416 RETURN
1417 ENDIF
1418 END DO
1419 ENDIF
1420 DO j = k+1, n
1421 m = augmatrix(j,k)/augmatrix(k,k)
1422 DO i = k, 2*n
1423 augmatrix(j,i) = augmatrix(j,i) - m*augmatrix(k,i)
1424 END DO
1425 END DO
1426 END DO
1427
1428c !TEST FOR INVERTIBILITY
1429 DO i = 1, n
1430 IF (augmatrix(i,i) == 0) THEN
1431 print*, "MATRIX IS NON - INVERTIBLE"
1432 inverse = 0
1433 errorflag = -1
1434 RETURN
1435 ENDIF
1436 END DO
1437
1438c !MAKE DIAGONAL ELEMENTS AS 1
1439 DO i = 1 , n
1440 m = augmatrix(i,i)
1441 DO j = i , (2 * n)
1442 augmatrix(i,j) = (augmatrix(i,j) / m)
1443 END DO
1444 END DO
1445
1446c !REDUCED RIGHT SIDE HALF OF AUGMENTED MATRIX TO IDENTITY MATRIX
1447 DO k = n-1, 1, -1
1448 DO i =1, k
1449 m = augmatrix(i,k+1)
1450 DO j = k, (2*n)
1451 augmatrix(i,j) = augmatrix(i,j) -augmatrix(k+1,j) * m
1452 END DO
1453 END DO
1454 END DO
1455
1456c !STORE ANSWER
1457 DO i =1, n
1458 DO j = 1, n
1459 inverse(i,j) = augmatrix(i,j+n)
1460 END DO
1461 END DO
1462 errorflag = 0
1463 RETURN

◆ prerbe3fr()

subroutine prerbe3fr ( integer, dimension(nrbe3l,*) irbe3,
integer n,
integer, dimension(3) jt,
integer, dimension(3) jr )

Definition at line 1512 of file hm_read_rbe3.F.

1513C-----------------------------------------------
1514C I m p l i c i t T y p e s
1515C-----------------------------------------------
1516#include "implicit_f.inc"
1517C-----------------------------------------------
1518C C o m m o n B l o c k s
1519C-----------------------------------------------
1520#include "param_c.inc"
1521C-----------------------------------------------
1522C D u m m y A r g u m e n t s
1523C-----------------------------------------------
1524 INTEGER IRBE3(NRBE3L,*),JT(3) ,JR(3),N
1525C REAL
1526C-----------------------------------------------
1527C L o c a l V a r i a b l e s
1528C-----------------------------------------------
1529 INTEGER J,IC,ICT,ICR
1530C======================================================================|
1531 ic=irbe3(4,n)
1532 ict=ic/512
1533 icr=(ic-512*(ict))/64
1534 DO j =1,3
1535 jt(j)=0
1536 jr(j)=0
1537 ENDDO
1538 SELECT CASE (ict)
1539 CASE(1)
1540 jt(3)=1
1541 CASE(2)
1542 jt(2)=1
1543 CASE(3)
1544 jt(2)=1
1545 jt(3)=1
1546 CASE(4)
1547 jt(1)=1
1548 CASE(5)
1549 jt(1)=1
1550 jt(3)=1
1551 CASE(6)
1552 jt(1)=1
1553 jt(2)=1
1554 CASE(7)
1555 jt(1)=1
1556 jt(2)=1
1557 jt(3)=1
1558 END SELECT
1559 SELECT CASE (icr)
1560 CASE(1)
1561 jr(3)=1
1562 CASE(2)
1563 jr(2)=1
1564 CASE(3)
1565 jr(2)=1
1566 jr(3)=1
1567 CASE(4)
1568 jr(1)=1
1569 CASE(5)
1570 jr(1)=1
1571 jr(3)=1
1572 CASE(6)
1573 jr(1)=1
1574 jr(2)=1
1575 CASE(7)
1576 jr(1)=1
1577 jr(2)=1
1578 jr(3)=1
1579 END SELECT
1580C---
1581 RETURN

◆ rbe3chk()

subroutine rbe3chk ( integer, dimension(*) inrbe3,
integer, dimension(*) ilrbe3,
integer ns,
xyz,
frbe3,
skew,
integer ng,
integer irot,
integer imodif,
wmin,
integer ipen,
integer ierr )

Definition at line 557 of file hm_read_rbe3.F.

560C-----------------------------------------------
561C I m p l i c i t T y p e s
562C-----------------------------------------------
563#include "implicit_f.inc"
564C-----------------------------------------------
565C C o m m o n B l o c k s
566C-----------------------------------------------
567#include "param_c.inc"
568C-----------------------------------------------
569C D u m m y A r g u m e n t s
570C-----------------------------------------------
571 INTEGER INRBE3(*),ILRBE3(*),NG, NS,IROT,IMODIF,IERR,IPEN
572C REAL
573 my_real
574 . xyz(3,*), frbe3(6,*), skew(lskew,*),wmin
575C-----------------------------------------------
576C L o c a l V a r i a b l e s
577C-----------------------------------------------
578 INTEGER I, J, K,KG,NSNGLR,IELSUB,KDIAG
579C REAL
580 my_real
581 * tw(3,ng), rw(3,ng),
582 * fufxlc(3,ng), fufylc(3,ng), fufzlc(3,ng),
583 * fumxlc(3,ng), fumylc(3,ng), fumzlc(3,ng),
584 * mxlc(3,ng), mylc(3,ng), mzlc(3,ng),
585 * fufx(3,ng), fufy(3,ng), fufz(3,ng),
586 * mufx(3,ng), mufy(3,ng), mufz(3,ng),
587 * fumx(3,ng), fumy(3,ng), fumz(3,ng),
588 * mx(3,ng), my(3,ng), mz(3,ng),
589 * mumx(3,ng), mumy(3,ng), mumz(3,ng),
590 * flocal(3,ng,6), mlocal(3,ng,6),
591 * fbasic(3,ng,6), mbasic(3,ng,6),
592 * fdstnl(3,ng,6), mdstnl(3,ng,6),
593 * fdstnb(3,ng,6), mdstnb(3,ng,6),el(3,3,ng)
594 my_real
595 * denfx, denfy, denfz, denmx, denmy, denmz,
596 * refpt(3), cgmx(3), cgmy(3), cgmz(3), averef,
597 * tfufx(3), tfufy(3), tfufz(3),
598 * tmufx(3), tmufy(3), tmufz(3),
599 * tfumx(3), tfumy(3), tfumz(3),
600 * tmumx(3), tmumy(3), tmumz(3),
601 * a(6,6), c(6,6), t(3,3),smin,smax,mmax,tmax,
602 * xbar(3),rn(3),gamma(9),wi(ng),gamma_max,rndotrn,det,arm
603C
604C INITIALIZATION
605C
606 CALL zero1(flocal,3*ng*6)
607 CALL zero1(mlocal,3*ng*6)
608 CALL zero1(fbasic,3*ng*6)
609 CALL zero1(mbasic,3*ng*6)
610 CALL zero1(fdstnl,3*ng*6)
611 CALL zero1(mdstnl,3*ng*6)
612 CALL zero1(fdstnb,3*ng*6)
613 CALL zero1(mdstnb,3*ng*6)
614 CALL zero1(a,36)
615 CALL zero1(c,36)
616 CALL zero1(cgmx,3)
617 CALL zero1(cgmy,3)
618 CALL zero1(cgmz,3)
619 ierr = 0
620 wmin =zero
621C----------debug use------------
622 kdiag = 0
623 refpt(1) = xyz(1,ns)
624 refpt(2) = xyz(2,ns)
625 refpt(3) = xyz(3,ns)
626 DO k = 1, ng
627 DO i = 1, 3
628 tw(i,k) = frbe3(i,k)
629 rw(i,k) = frbe3(i+3,k)
630 ENDDO
631 ENDDO
632C
633C ERROR OUT IF RBE3 ELEMENT HAS TWO INDEPENDENT NODES WITH
634C NO ROTATIONAL WEIGHTS SET (THIS MEANS THE ELEMENT CANNOT
635C SUPPORT A MOMENT ALONG ITS AXIS)
636C
637 IF (ng == 2.AND.irot==0) THEN
638 ierr = 322
639 GOTO 999
640 ENDIF
641C
642C CALCULATE DIRECTION COSINES OF LOCAL COORDINATE SYSTEMS, IF ANY
643C
644 DO k = 1, ng
645 ielsub = ilrbe3(k)
646 IF (ielsub > 0) THEN
647 DO i = 1, 3
648 el(i,1,k) = skew(i,ielsub)
649 el(i,2,k) = skew(i+3,ielsub)
650 el(i,3,k) = skew(i+6,ielsub)
651 ENDDO
652 ENDIF
653 ENDDO
654C
655C DENOMINATORS FOR DISTRIBUTING FORCES (DENFX, DENFY AND DENFZ)
656C
657 denfx = zero
658 denfy = zero
659 denfz = zero
660 averef = zero
661C
662 DO 70 k = 1, ng
663 kg = inrbe3(k)
664 ielsub = ilrbe3(k)
665 IF (ielsub > 0) THEN
666C
667C IF GRID POINT HAS A LOCAL COORDINATE SYSTEM
668C
669 DO 60 i = 1, 3
670 denfx = denfx + tw(i,k)*el(i,1,k)**2
671 denfy = denfy + tw(i,k)*el(i,2,k)**2
672 denfz = denfz + tw(i,k)*el(i,3,k)**2
673 60 CONTINUE
674 ELSE
675 denfx = denfx + tw(1,k)
676 denfy = denfy + tw(2,k)
677 denfz = denfz + tw(3,k)
678 END IF
679C
680 averef = averef + sqrt( (xyz(1,kg) - refpt(1))**2 +
681 * (xyz(2,kg) - refpt(2))**2 +
682 * (xyz(3,kg) - refpt(3))**2 )
683 70 CONTINUE
684C
685 IF (abs(denfx) <= em20) THEN
686 ierr = 326
687 ENDIF
688C
689 IF (abs(denfy) <= em20) THEN
690 ierr = 327
691 ENDIF
692C
693 IF (abs(denfz) <= em20) THEN
694 ierr = 328
695 ENDIF
696 IF (ierr > 0) GOTO 999
697 averef = averef/ng
698 IF (averef == zero) averef = 1.0d0
699C--- IMODIF=4: normalized TW;
700 IF (imodif==4.OR.ipen>0) THEN
701 DO k = 1, ng
702 frbe3(1,k) = frbe3(1,k)/denfx
703 frbe3(2,k) = frbe3(2,k)/denfy
704 frbe3(3,k) = frbe3(3,k)/denfz
705 frbe3(4,k) = frbe3(4,k)/denfx
706 frbe3(5,k) = frbe3(5,k)/denfy
707 frbe3(6,k) = frbe3(6,k)/denfz
708 ENDDO
709 END IF
710 IF (ipen > 0) THEN
711 xbar(1:3) = zero
712 DO k = 1, ng
713 kg = inrbe3(k)
714 wi(k) = frbe3(1,k)
715 xbar(1:3) = xbar(1:3) + wi(k)*xyz(1:3,kg)
716 ENDDO
717! keep kinematic if very small arm
718 rn(1:3) = refpt(1:3)-xbar(1:3)
719 arm = rn(1)*rn(1)+rn(2)*rn(2)+rn(3)*rn(3)
720 rndotrn = zero
721 DO k = 1, ng
722 kg = inrbe3(k)
723 rn(1:3) = xyz(1:3,kg)-xbar(1:3)
724 rndotrn = max(rndotrn,rn(1)*rn(1)+rn(2)*rn(2)+rn(3)*rn(3))
725 END DO
726 IF (arm/rndotrn < em06) ipen =-2
727C--- CHECK for penalty colinear case w/ irot=0
728 IF (irot==0) THEN
729 gamma(1:9) = zero
730 DO k = 1, ng
731 kg = inrbe3(k)
732 rn(1:3) = xyz(1:3,kg)-xbar(1:3)
733 rndotrn = rn(1)*rn(1)+rn(2)*rn(2)+rn(3)*rn(3)
734!
735 gamma(1) = gamma(1)+wi(k)*(rndotrn-rn(1)*rn(1))
736 gamma(2) = gamma(2)+wi(k)*( -rn(2)*rn(1))
737 gamma(3) = gamma(3)+wi(k)*( -rn(3)*rn(1))
738 gamma(4) = gamma(4)+wi(k)*( -rn(1)*rn(2))
739 gamma(5) = gamma(5)+wi(k)*(rndotrn-rn(2)*rn(2))
740 gamma(6) = gamma(6)+wi(k)*( -rn(3)*rn(2))
741 gamma(7) = gamma(7)+wi(k)*( -rn(1)*rn(3))
742 gamma(8) = gamma(8)+wi(k)*( -rn(2)*rn(3))
743 gamma(9) = gamma(9)+wi(k)*(rndotrn-rn(3)*rn(3))
744 ENDDO
745 det = (gamma(1)*(gamma(5)*gamma(9)-gamma(6)*gamma(8))-
746 * gamma(2)*(gamma(4)*gamma(9)-gamma(6)*gamma(7))+
747 * gamma(3)*(gamma(4)*gamma(8)-gamma(5)*gamma(7)))
748!
749 gamma_max = max(em20,gamma(1),gamma(5),gamma(9))
750 IF(abs(det/(gamma_max*gamma_max*gamma_max)) < em6) ierr = 400
751 END IF
752 END IF
753 IF (ierr > 0) GOTO 999
754!
755C
756C CALCULATE 3 CENTERS OF GRAVITY (CGMX, CGMY AND CGMZ) AND
757C DENOMINATORS FOR DISTRIBUTING MOMENTS (DENMX, DENMY AND DENMZ)
758C
759 DO 40 k = 1, ng
760 kg = inrbe3(k)
761 ielsub = ilrbe3(k)
762 IF (ielsub > 0) THEN
763C
764C IF THERE IS A LOCAL COORDINATE SYSTEM AT THE GRID POINT
765C
766 DO 10 i = 1, 3
767 cgmx(2) = cgmx(2) + tw(i,k)*el(i,3,k)**2*xyz(2,kg)
768 cgmx(3) = cgmx(3) + tw(i,k)*el(i,2,k)**2*xyz(3,kg)
769 10 CONTINUE
770C
771 DO 20 i = 1, 3
772 cgmy(3) = cgmy(3) + tw(i,k)*el(i,1,k)**2*xyz(3,kg)
773 cgmy(1) = cgmy(1) + tw(i,k)*el(i,3,k)**2*xyz(1,kg)
774 20 CONTINUE
775C
776 DO 30 i = 1, 3
777 cgmz(1) = cgmz(1) + tw(i,k)*el(i,2,k)**2*xyz(1,kg)
778 cgmz(2) = cgmz(2) + tw(i,k)*el(i,1,k)**2*xyz(2,kg)
779 30 CONTINUE
780C
781 ELSE
782 cgmx(2) = cgmx(2) + tw(3,k)*xyz(2,kg)
783 cgmx(3) = cgmx(3) + tw(2,k)*xyz(3,kg)
784C
785 cgmy(3) = cgmy(3) + tw(1,k)*xyz(3,kg)
786 cgmy(1) = cgmy(1) + tw(3,k)*xyz(1,kg)
787C
788 cgmz(1) = cgmz(1) + tw(2,k)*xyz(1,kg)
789 cgmz(2) = cgmz(2) + tw(1,k)*xyz(2,kg)
790 END IF
791 40 CONTINUE
792 cgmx(2) = cgmx(2)/denfz
793 cgmx(3) = cgmx(3)/denfy
794C
795 cgmy(3) = cgmy(3)/denfx
796 cgmy(1) = cgmy(1)/denfz
797C
798 cgmz(1) = cgmz(1)/denfy
799 cgmz(2) = cgmz(2)/denfx
800C
801 denmx = zero
802 denmy = zero
803 denmz = zero
804C
805 DO 90 k = 1, ng
806 kg = inrbe3(k)
807 ielsub = ilrbe3(k)
808C
809C NOTE: AS IMPLEMENTED IN NASTRAN 70.7, WE SCALE THE ROTATIONAL
810C WEIGHTS WITH THE SQUARE OF THE AVERAGE DISTANCE OF THE
811C INDEPENDENT GRID POINTS FROM THE REFERENCE POINT TO
812C RENDER THE RBE3 CALCULATIONS UNIT INDEPENDENT
813C
814 IF (ielsub > 0) THEN
815C
816C IF GRID POINT HAS A LOCAL COORDINATE SYSTEM
817C
818 DO 80 i = 1, 3
819 denmx = denmx + rw(i,k)*el(i,1,k)**2*averef**2 +
820 * tw(i,k)*( el(i,3,k)*(xyz(2,kg) - cgmx(2)) -
821 * el(i,2,k)*(xyz(3,kg) - cgmx(3))
822 * ) **2
823 denmy = denmy + rw(i,k)*el(i,2,k)**2*averef**2 +
824 * tw(i,k)*( el(i,1,k)*(xyz(3,kg) - cgmy(3)) -
825 * el(i,3,k)*(xyz(1,kg) - cgmy(1))
826 * ) **2
827 denmz = denmz + rw(i,k)*el(i,3,k)**2*averef**2 +
828 * tw(i,k)*( el(i,2,k)*(xyz(1,kg) - cgmz(1)) -
829 * el(i,1,k)*(xyz(2,kg) - cgmz(2))
830 * ) **2
831 80 CONTINUE
832 ELSE
833 denmx = denmx + rw(1,k)*averef**2 +
834 * tw(2,k)*(xyz(3,kg) - cgmx(3))**2 +
835 * tw(3,k)*(xyz(2,kg) - cgmx(2))**2
836 denmy = denmy + rw(2,k)*averef**2 +
837 * tw(1,k)*(xyz(3,kg) - cgmy(3))**2 +
838 * tw(3,k)*(xyz(1,kg) - cgmy(1))**2
839 denmz = denmz + rw(3,k)*averef**2 +
840 * tw(2,k)*(xyz(1,kg) - cgmz(1))**2 +
841 * tw(1,k)*(xyz(2,kg) - cgmz(2))**2
842 END IF
843 90 CONTINUE
844C
845C PERFORM SOME CHECKS ON WEIGHTS, TO MAKE SURE THAT THE RBE3
846C ELEMENT HAS NO UNCONSTRAINED DEGREES OF FREEDOM
847C
848C
849 IF (abs(denmx) <= em20) THEN
850 ierr = 329
851 ENDIF
852C
853 IF (abs(denmy) <= em20) THEN
854 ierr = 330
855 ENDIF
856C
857 IF (abs(denmz) <= em20) THEN
858 ierr = 331
859 ENDIF
860C
861 smin = min(abs(denmx),abs(denmy),abs(denmz))
862 smax = max(abs(denmx),abs(denmy),abs(denmz))
863C
864 IF (ierr > 0) GOTO 999
865C
866 IF (irot==0 .AND.(smax/smin)>thirty) ierr = -100
867C CALCULATE 3 FORCE DISTRIBUTIONS THAT CREATE NET X, Y AND Z FORCES
868C OF 1 (BESIDES OTHER NONZERO FORCES/MOMENTS IN ALL THE DIRECTIONS)
869C
870 CALL rbe3uf(inrbe3,ilrbe3,el,tw,xyz,refpt,
871 * fufxlc,fufylc,fufzlc,fufx,fufy,fufz,mufx,mufy,mufz,
872 * tfufx,tfufy,tfufz,tmufx,tmufy,tmufz,
873 * denfx,denfy,denfz,ng)
874C
875C CALCULATE 3 MOMENT/FORCE DISTRIBUTIONS THAT CREATE NET X, Y AND Z
876C MOMENTS OF 1 (BESIDES OTHER NONZERO FORCES/MOMENTS IN ALL THE
877C DIRECTIONS) AT CGMX, CGMY AND CGMZ RESPECTIVELY
878C
879 CALL rbe3um(inrbe3,ilrbe3,el,tw,rw,xyz,refpt,cgmx,cgmy,cgmz,
880 * fumxlc,fumylc,fumzlc,mxlc,mylc,mzlc,
881 * fumx,fumy,fumz,mx,my,mz,mumx,mumy,mumz,
882 * tfumx,tfumy,tfumz,tmumx,tmumy,tmumz,
883 * averef,denmx,denmy,denmz,ng,irot )
884C
885C DETERMINE COMBINATORY COEFFICIENTS FOR THESE 6 DISTRIBUTIONS
886C (6 COEFFICIENTS FOR EACH OF 6 CASES)
887C
888C CASE 1 - RESULTANT OF THE LINEAR COMBINATION OF THESE FORCE/MOMENT
889C DISTRIBUTIONS IS A UNIT X-FORCE AT REFERENCE POINT
890C CASE 2 - RESULTANT OF THE LINEAR COMBINATION OF THESE FORCE/MOMENT
891C DISTRIBUTIONS IS A UNIT Y-FORCE AT REFERENCE POINT
892C CASE 3 - RESULTANT OF THE LINEAR COMBINATION OF THESE FORCE/MOMENT
893C DISTRIBUTIONS IS A UNIT Z-FORCE AT REFERENCE POINT
894C CASE 4 - RESULTANT OF THE LINEAR COMBINATION OF THESE FORCE/MOMENT
895C DISTRIBUTIONS IS A UNIT X-MOMENT AT REFERENCE POINT
896C CASE 5 - RESULTANT OF THE LINEAR COMBINATION OF THESE FORCE/MOMENT
897C DISTRIBUTIONS IS A UNIT Y-MOMENT AT REFERENCE POINT
898C CASE 6 - RESULTANT OF THE LINEAR COMBINATION OF THESE FORCE/MOMENT
899C DISTRIBUTIONS IS A UNIT Z-MOMENT AT REFERENCE POINT
900C
901C IN ORDER TO DETERMINE THESE COEFFICIENTS, FIRST SET UP A 6X6
902C MATRIX. THE 6 COLUMNS OF THE INVERSE OF THIS MATRIX ARE THE
903C DESIRED 6 SETS OF COEFFICIENTS.
904C
905 DO 120 i = 1, 3
906 k = i + 3
907 a(i,1) = tfufx(i)
908 a(k,1) = tmufx(i)
909 a(i,2) = tfufy(i)
910 a(k,2) = tmufy(i)
911 a(i,3) = tfufz(i)
912 a(k,3) = tmufz(i)
913 a(i,4) = tfumx(i)
914 a(k,4) = tmumx(i)
915 a(i,5) = tfumy(i)
916 a(k,5) = tmumy(i)
917 a(i,6) = tfumz(i)
918 a(k,6) = tmumz(i)
919 120 CONTINUE
920C
921C INVERT THE 6X6 MATRIX
922C
923 nsnglr = 0
924 CALL invert(a,c,6,nsnglr)
925 IF (nsnglr /= 0) THEN
926 ierr = 332
927 GOTO 999
928 ENDIF
929 IF (kdiag >= 1) THEN
930 CALL wrrinf('C(i,1)=',c(1,1),3)
931 CALL wrrinf('C(i,2)=',c(1,2),3)
932 CALL wrrinf('C(i,3)=',c(1,3),3)
933 ENDIF
934 IF (kdiag==0.AND.ierr==0) RETURN
935C
936C LINEARLY COMBINE FORCE/MOMENT DISTRIBUTIONS FOR THE 6 CASES
937C
938 DO 170 j = 1, 6
939 DO 160 k = 1, ng
940 DO 150 i = 1, 3
941 flocal(i,k,j) = c(1,j)*fufxlc(i,k) + c(2,j)*fufylc(i,k) +
942 * c(3,j)*fufzlc(i,k) + c(4,j)*fumxlc(i,k) +
943 * c(5,j)*fumylc(i,k) + c(6,j)*fumzlc(i,k)
944 mlocal(i,k,j) = c(4,j)*mxlc(i,k) + c(5,j)*mylc(i,k) +
945 * c(6,j)*mzlc(i,k)
946 fbasic(i,k,j) = c(1,j)*fufx(i,k) + c(2,j)*fufy(i,k) +
947 * c(3,j)*fufz(i,k) + c(4,j)*fumx(i,k) +
948 * c(5,j)*fumy(i,k) + c(6,j)*fumz(i,k)
949 mbasic(i,k,j) = c(4,j)*mx(i,k) + c(5,j)*my(i,k) +
950 * c(6,j)*mz(i,k)
951 150 CONTINUE
952 160 CONTINUE
953 170 CONTINUE
954C
955C
956C NO LOCAL COORDINATE SYSTEM AT REFERENCE POINT
957C
958 DO 270 j = 1, 6
959 DO 260 k = 1, ng
960 DO 250 i = 1, 3
961 fdstnl(i,k,j) = flocal(i,k,j)
962 mdstnl(i,k,j) = mlocal(i,k,j)
963 fdstnb(i,k,j) = fbasic(i,k,j)
964 mdstnb(i,k,j) = mbasic(i,k,j)
965 250 CONTINUE
966 260 CONTINUE
967 270 CONTINUE
968C--------------special case with Imodif
969 IF (ierr==-100) THEN
970 mmax=zero
971 DO j = 4, 6
972 DO k = 1, ng
973 DO i = 1, 3
974 IF (mmax<abs(fdstnb(i,k,j))) mmax = abs(fdstnb(i,k,j))
975 END DO
976 END DO
977 END DO
978 IF (mmax<=one) THEN
979 ierr=0
980 ELSE
981 tmax=zero
982 IF (imodif/=2) THEN
983 DO k = 1, ng
984 DO i = 1, 3
985 IF (tmax<tw(i,k)) tmax=tw(i,k)
986 ENDDO
987 ENDDO
988 ENDIF
989 wmin=tmax*em04
990 DO k = 1, ng
991 frbe3(1,k) = max(wmin,frbe3(1,k))
992 frbe3(2,k) = max(wmin,frbe3(2,k))
993 frbe3(3,k) = max(wmin,frbe3(3,k))
994 ENDDO
995 ENDIF
996 END IF
997C
998 999 CONTINUE
999C
1000C DIAGNOSTIC INFORMATION
1001C
1002 IF (kdiag >= 2) THEN
1003c CALL WRRINF('REF_POINT',REFPT,3)
1004c CALL WRRINF('INDPT_GRDS',INRBE3,NG)
1005c IF (SMAX/SMIN > THIRTY) print *,'SMAX/SMIN=',SMAX/SMIN
1006 CALL wrrinf('TRAN_WGHTS',tw,3*ng)
1007 CALL wrrinf('ROT_WGHTS',rw,3*ng)
1008 CALL wrrinf('CGMX',cgmx,3)
1009 CALL wrrinf('CGMY',cgmy,3)
1010 CALL wrrinf('CGMZ',cgmz,3)
1011 CALL wrrinf('DENFX',denfx,1)
1012 CALL wrrinf('DENFY',denfy,1)
1013 CALL wrrinf('DENFZ',denfz,1)
1014 CALL wrrinf('DENMX',denmx,1)
1015 CALL wrrinf('DENMY',denmy,1)
1016 CALL wrrinf('DENMZ',denmz,1)
1017 CALL wrrinf('AVEREF',averef,1)
1018C
1019 IF (kdiag == 9.or.ierr/=0) THEN
1020 CALL wrrinf('FDSTNB_ULFX@REF',fdstnb(1,1,1),3*ng)
1021 CALL wrrinf('FDSTNB_ULFY@REF',fdstnb(1,1,2),3*ng)
1022 CALL wrrinf('FDSTNB_ULFZ@REF',fdstnb(1,1,3),3*ng)
1023 CALL wrrinf('FDSTNB_ULMX@REF',fdstnb(1,1,4),3*ng)
1024 CALL wrrinf('FDSTNB_ULMY@REF',fdstnb(1,1,5),3*ng)
1025 CALL wrrinf('FDSTNB_ULMZ@REF',fdstnb(1,1,6),3*ng)
1026 CALL wrrinf('MDSTNB_ULFX@REF',mdstnb(1,1,1),3*ng)
1027 CALL wrrinf('MDSTNB_ULFY@REF',mdstnb(1,1,2),3*ng)
1028 CALL wrrinf('MDSTNB_ULFZ@REF',mdstnb(1,1,3),3*ng)
1029 CALL wrrinf('MDSTNB_ULMX@REF',mdstnb(1,1,4),3*ng)
1030 CALL wrrinf('MDSTNB_ULMY@REF',mdstnb(1,1,5),3*ng)
1031 CALL wrrinf('MDSTNB_ULMZ@REF',mdstnb(1,1,6),3*ng)
1032 ENDIF
1033 IF (kdiag >= 30) THEN
1034 CALL wrrinf('FDSTNL_ULFX@REF',fdstnl(1,1,1),3*ng)
1035 CALL wrrinf('FDSTNL_ULFY@REF',fdstnl(1,1,2),3*ng)
1036 CALL wrrinf('FDSTNL_ULFZ@REF',fdstnl(1,1,3),3*ng)
1037 CALL wrrinf('FDSTNL_ULMX@REF',fdstnl(1,1,4),3*ng)
1038 CALL wrrinf('FDSTNL_ULMY@REF',fdstnl(1,1,5),3*ng)
1039 CALL wrrinf('FDSTNL_ULMZ@REF',fdstnl(1,1,6),3*ng)
1040 CALL wrrinf('MDSTNL_ULFX@REF',mdstnl(1,1,1),3*ng)
1041 CALL wrrinf('MDSTNL_ULFY@REF',mdstnl(1,1,2),3*ng)
1042 CALL wrrinf('MDSTNL_ULFZ@REF',mdstnl(1,1,3),3*ng)
1043 CALL wrrinf('MDSTNL_ULMX@REF',mdstnl(1,1,4),3*ng)
1044 CALL wrrinf('MDSTNL_ULMY@REF',mdstnl(1,1,5),3*ng)
1045 CALL wrrinf('MDSTNL_ULMZ@REF',mdstnl(1,1,6),3*ng)
1046C
1047C
1048 ENDIF
1049 ENDIF
1050C
1051 RETURN
subroutine invert(matrix, inverse, n, errorflag)
subroutine zero1(r, n)
subroutine wrrinf(title, r, n)
subroutine rbe3um(inrbe3, ilrbe3, el, tw, rw, xyz, refpt, cgmx, cgmy, cgmz, fumxlc, fumylc, fumzlc, mxlc, mylc, mzlc, fumx, fumy, fumz, mx, my, mz, mumx, mumy, mumz, tfumx, tfumy, tfumz, tmumx, tmumy, tmumz, averef, denmx, denmy, denmz, ng, irot)
subroutine rbe3uf(inrbe3, ilrbe3, el, tw, xyz, refpt, fufxlc, fufylc, fufzlc, fufx, fufy, fufz, mufx, mufy, mufz, tfufx, tfufy, tfufz, tmufx, tmufy, tmufz, denfx, denfy, denfz, ng)
#define max(a, b)
Definition macros.h:21

◆ rbe3uf()

subroutine rbe3uf ( integer, dimension(ng) inrbe3,
integer, dimension(ng) ilrbe3,
el,
tw,
xyz,
refpt,
fufxlc,
fufylc,
fufzlc,
fufx,
fufy,
fufz,
mufx,
mufy,
mufz,
tfufx,
tfufy,
tfufz,
tmufx,
tmufy,
tmufz,
denfx,
denfy,
denfz,
integer ng )

Definition at line 1061 of file hm_read_rbe3.F.

1066C-----------------------------------------------
1067C I m p l i c i t T y p e s
1068C-----------------------------------------------
1069#include "implicit_f.inc"
1070 INTEGER NG
1071 INTEGER INRBE3(NG), ILRBE3(NG)
1072 my_real
1073 * el(3,3,*),tw(3,ng), xyz(3,*), refpt(3),
1074 * fufxlc(3,ng), fufylc(3,ng), fufzlc(3,ng),
1075 * fufx(3,ng), fufy(3,ng), fufz(3,ng),
1076 * mufx(3,ng), mufy(3,ng), mufz(3,ng),
1077 * tfufx(3), tfufy(3), tfufz(3),
1078 * tmufx(3), tmufy(3), tmufz(3)
1079 my_real
1080 * denfx, denfy, denfz,xarm, yarm, zarm
1081 INTEGER I, J, K, KG, IELSUB
1082C
1083C INITIALIZE FORCE AND MOMENT DISTRIBUTIONS TO ZERO
1084C
1085 CALL zero1(fufx,3*ng)
1086 CALL zero1(fufy,3*ng)
1087 CALL zero1(fufz,3*ng)
1088 CALL zero1(tfufx,3)
1089 CALL zero1(tfufy,3)
1090 CALL zero1(tfufz,3)
1091 CALL zero1(tmufx,3)
1092 CALL zero1(tmufy,3)
1093 CALL zero1(tmufz,3)
1094C
1095C FORCE DISTRIBUTIONS AT RBE3 GRID POINTS CORRESPONDING TO UNIT
1096C APPLIED FORCES AT RBE3 REFERENCE POINT ALONG (BASIC COORDINATE)
1097C X, Y AND Z DIRECTIONS
1098C
1099 DO 50 k = 1, ng
1100 kg = inrbe3(k)
1101 ielsub = ilrbe3(k)
1102 IF (ielsub > 0) THEN
1103C
1104C FORCES AT GRID POINT ALONG GRID POINT'S LOCAL (OUTPUT)
1105C COORDINATE AXES
1106C
1107 DO 10 i = 1, 3
1108 fufxlc(i,k) = tw(i,k)*el(i,1,k)/denfx
1109 fufylc(i,k) = tw(i,k)*el(i,2,k)/denfy
1110 fufzlc(i,k) = tw(i,k)*el(i,3,k)/denfz
1111 10 CONTINUE
1112C
1113C FORCES AT GRID POINT ALONG BASIC COORDINATE AXES
1114C
1115 DO 30 i = 1, 3
1116 DO 20 j = 1, 3
1117 fufx(j,k) = fufx(j,k) + fufxlc(i,k)*el(i,j,k)
1118 fufy(j,k) = fufy(j,k) + fufylc(i,k)*el(i,j,k)
1119 fufz(j,k) = fufz(j,k) + fufzlc(i,k)*el(i,j,k)
1120 20 CONTINUE
1121 30 CONTINUE
1122C
1123 ELSE
1124 fufxlc(1,k) = tw(1,k)/denfx
1125 fufylc(2,k) = tw(2,k)/denfy
1126 fufzlc(3,k) = tw(3,k)/denfz
1127 fufx(1,k) = fufxlc(1,k)
1128 fufy(2,k) = fufylc(2,k)
1129 fufz(3,k) = fufzlc(3,k)
1130 ENDIF
1131C
1132C MOMENTS AT REFERENCE POINT DUE TO THESE FORCE DISTRIBUTIONS
1133C
1134 xarm = xyz(1,kg) - refpt(1)
1135 yarm = xyz(2,kg) - refpt(2)
1136 zarm = xyz(3,kg) - refpt(3)
1137C
1138C MOMENTS AT REFERENCE POINT DUE TO FUFX
1139C
1140 mufx(1,k) = yarm*fufx(3,k) - zarm*fufx(2,k)
1141 mufx(2,k) = zarm*fufx(1,k) - xarm*fufx(3,k)
1142 mufx(3,k) = xarm*fufx(2,k) - yarm*fufx(1,k)
1143C
1144C MOMENTS AT REFERENCE POINT DUE TO FUFY
1145C
1146 mufy(1,k) = yarm*fufy(3,k) - zarm*fufy(2,k)
1147 mufy(2,k) = zarm*fufy(1,k) - xarm*fufy(3,k)
1148 mufy(3,k) = xarm*fufy(2,k) - yarm*fufy(1,k)
1149C
1150C MOMENTS AT REFERENCE POINT DUE TO FUFZ
1151C
1152 mufz(1,k) = yarm*fufz(3,k) - zarm*fufz(2,k)
1153 mufz(2,k) = zarm*fufz(1,k) - xarm*fufz(3,k)
1154 mufz(3,k) = xarm*fufz(2,k) - yarm*fufz(1,k)
1155C
1156C TOTAL FORCES AND MOMENTS
1157C
1158 DO 40 j = 1, 3
1159 tfufx(j) = tfufx(j) + fufx(j,k)
1160 tfufy(j) = tfufy(j) + fufy(j,k)
1161 tfufz(j) = tfufz(j) + fufz(j,k)
1162 tmufx(j) = tmufx(j) + mufx(j,k)
1163 tmufy(j) = tmufy(j) + mufy(j,k)
1164 tmufz(j) = tmufz(j) + mufz(j,k)
1165 40 CONTINUE
1166C
1167 50 CONTINUE
1168C
1169 RETURN

◆ rbe3um()

subroutine rbe3um ( integer, dimension(ng) inrbe3,
integer, dimension(ng) ilrbe3,
el,
tw,
rw,
xyz,
refpt,
cgmx,
cgmy,
cgmz,
fumxlc,
fumylc,
fumzlc,
mxlc,
mylc,
mzlc,
fumx,
fumy,
fumz,
mx,
my,
mz,
mumx,
mumy,
mumz,
tfumx,
tfumy,
tfumz,
tmumx,
tmumy,
tmumz,
averef,
denmx,
denmy,
denmz,
integer ng,
integer irot )

Definition at line 1180 of file hm_read_rbe3.F.

1185C-----------------------------------------------
1186C I m p l i c i t T y p e s
1187C-----------------------------------------------
1188#include "implicit_f.inc"
1189 INTEGER NG,IROT
1190 INTEGER INRBE3(NG), ILRBE3(NG)
1191 my_real
1192 * el(3,3,*),tw(3,ng), rw(3,ng), xyz(3,*),
1193 * refpt(3), cgmx(3), cgmy(3), cgmz(3),
1194 * fumxlc(3,ng), fumylc(3,ng), fumzlc(3,ng),
1195 * mxlc(3,ng), mylc(3,ng), mzlc(3,ng),
1196 * fumx(3,ng), fumy(3,ng), fumz(3,ng),
1197 * mx(3,ng), my(3,ng), mz(3,ng),
1198 * mumx(3,ng), mumy(3,ng), mumz(3,ng),
1199 * tfumx(3), tfumy(3), tfumz(3),
1200 * tmumx(3), tmumy(3), tmumz(3)
1201 my_real
1202 * averef, denmx, denmy, denmz,xarm, yarm, zarm
1203 INTEGER I, J, K, KG, IELSUB
1204C
1205C INITIALIZE FORCE AND MOMENT DISTRIBUTIONS TO ZERO
1206C
1207 CALL zero1(fumx,3*ng)
1208 CALL zero1(fumy,3*ng)
1209 CALL zero1(fumz,3*ng)
1210 CALL zero1(mx,3*ng)
1211 CALL zero1(my,3*ng)
1212 CALL zero1(mz,3*ng)
1213 CALL zero1(tfumx,3)
1214 CALL zero1(tfumy,3)
1215 CALL zero1(tfumz,3)
1216 CALL zero1(tmumx,3)
1217 CALL zero1(tmumy,3)
1218 CALL zero1(tmumz,3)
1219C
1220C FORCE AND MOMENT DISTRIBUTIONS AT RBE3 GRID POINTS CORRESPONDING
1221C TO UNIT APPLIED MOMENTS AT RBE3 REFERENCE POINT ALONG (BASIC
1222C COORDINATE) X, Y AND Z DIRECTIONS
1223C
1224 DO 50 k = 1, ng
1225 kg = inrbe3(k)
1226 ielsub = ilrbe3(k)
1227 IF (ielsub > 0) THEN
1228C
1229C FORCES AT GRID POINT ALONG GRID POINT'S LOCAL
1230C (OUTPUT) COORDINATE AXES
1231C
1232 DO 10 i = 1, 3
1233 fumxlc(i,k) = tw(i,k)*
1234 * ( el(i,3,k)*(xyz(2,kg) - cgmx(2)) -
1235 * el(i,2,k)*(xyz(3,kg) - cgmx(3))
1236 * )/denmx
1237 fumylc(i,k) = tw(i,k)*
1238 * ( el(i,1,k)*(xyz(3,kg) - cgmy(3)) -
1239 * el(i,3,k)*(xyz(1,kg) - cgmy(1))
1240 * )/denmy
1241 fumzlc(i,k) = tw(i,k)*
1242 * ( el(i,2,k)*(xyz(1,kg) - cgmz(1)) -
1243 * el(i,1,k)*(xyz(2,kg) - cgmz(2))
1244 * )/denmz
1245 10 CONTINUE
1246C
1247C FORCES AND MOMENTS AT GRID POINT ALONG BASIC COORDINATE AXES
1248C
1249 DO 30 i = 1, 3
1250 DO 20 j = 1, 3
1251 fumx(j,k) = fumx(j,k) + fumxlc(i,k)*el(i,j,k)
1252 fumy(j,k) = fumy(j,k) + fumylc(i,k)*el(i,j,k)
1253 fumz(j,k) = fumz(j,k) + fumzlc(i,k)*el(i,j,k)
1254 20 CONTINUE
1255 30 CONTINUE
1256C
1257 ELSE
1258 fumxlc(2,k) = -tw(2,k)*(xyz(3,kg) - cgmx(3))/denmx
1259 fumxlc(3,k) = tw(3,k)*(xyz(2,kg) - cgmx(2))/denmx
1260 fumylc(1,k) = tw(1,k)*(xyz(3,kg) - cgmy(3))/denmy
1261 fumylc(3,k) = -tw(3,k)*(xyz(1,kg) - cgmy(1))/denmy
1262 fumzlc(1,k) = -tw(1,k)*(xyz(2,kg) - cgmz(2))/denmz
1263 fumzlc(2,k) = tw(2,k)*(xyz(1,kg) - cgmz(1))/denmz
1264C
1265 fumx(2,k) = fumxlc(2,k)
1266 fumx(3,k) = fumxlc(3,k)
1267 fumy(1,k) = fumylc(1,k)
1268 fumy(3,k) = fumylc(3,k)
1269 fumz(1,k) = fumzlc(1,k)
1270 fumz(2,k) = fumzlc(2,k)
1271 ENDIF
1272C
1273C MOMENTS AT REFERENCE POINT DUE TO FUMX
1274C
1275 xarm = xyz(1,kg) - refpt(1)
1276 yarm = xyz(2,kg) - refpt(2)
1277 zarm = xyz(3,kg) - refpt(3)
1278C
1279 mumx(1,k) = yarm*fumx(3,k) - zarm*fumx(2,k)
1280 mumx(2,k) = zarm*fumx(1,k) - xarm*fumx(3,k)
1281 mumx(3,k) = xarm*fumx(2,k) - yarm*fumx(1,k)
1282C
1283C MOMENTS AT REFERENCE POINT DUE TO FUMY
1284C
1285 mumy(1,k) = yarm*fumy(3,k) - zarm*fumy(2,k)
1286 mumy(2,k) = zarm*fumy(1,k) - xarm*fumy(3,k)
1287 mumy(3,k) = xarm*fumy(2,k) - yarm*fumy(1,k)
1288C
1289C MOMENTS AT REFERENCE POINT DUE TO FUMZ
1290C
1291 mumz(1,k) = yarm*fumz(3,k) - zarm*fumz(2,k)
1292 mumz(2,k) = zarm*fumz(1,k) - xarm*fumz(3,k)
1293 mumz(3,k) = xarm*fumz(2,k) - yarm*fumz(1,k)
1294C
1295 50 CONTINUE
1296C
1297 IF (irot>0) THEN
1298 DO k = 1, ng
1299 kg = inrbe3(k)
1300 ielsub = ilrbe3(k)
1301 IF (ielsub > 0) THEN
1302C
1303C MOMENTS AT GRID POINT ALONG GRID POINT'S LOCAL
1304C (OUTPUT) COORDINATE AXES
1305C
1306 DO i = 1, 3
1307 mxlc(i,k) = averef**2*rw(i,k)*el(i,1,k)/denmx
1308 mylc(i,k) = averef**2*rw(i,k)*el(i,2,k)/denmy
1309 mzlc(i,k) = averef**2*rw(i,k)*el(i,3,k)/denmz
1310 END DO
1311C
1312C MOMENTS AT GRID POINT ALONG BASIC COORDINATE AXES
1313C
1314 DO i = 1, 3
1315 DO j = 1, 3
1316 mx(j,k) = mx(j,k) + mxlc(i,k)*el(i,j,k)
1317 my(j,k) = my(j,k) + mylc(i,k)*el(i,j,k)
1318 mz(j,k) = mz(j,k) + mzlc(i,k)*el(i,j,k)
1319 END DO
1320 END DO
1321C
1322 ELSE
1323 mxlc(1,k) = averef**2*rw(1,k)/denmx
1324 mylc(2,k) = averef**2*rw(2,k)/denmy
1325 mzlc(3,k) = averef**2*rw(3,k)/denmz
1326C
1327 mx(1,k) = mxlc(1,k)
1328 my(2,k) = mylc(2,k)
1329 mz(3,k) = mzlc(3,k)
1330 ENDIF
1331C
1332 DO j = 1, 3
1333 mumx(j,k) = mumx(j,k) + mx(j,k)
1334 mumy(j,k) = mumy(j,k) + my(j,k)
1335 mumz(j,k) = mumz(j,k) + mz(j,k)
1336 END DO
1337 END DO
1338 END IF
1339C
1340C
1341C TOTAL FORCES AND MOMENTS
1342C
1343C
1344 DO k = 1, ng
1345 DO j = 1, 3
1346 tfumx(j) = tfumx(j) + fumx(j,k)
1347 tfumy(j) = tfumy(j) + fumy(j,k)
1348 tfumz(j) = tfumz(j) + fumz(j,k)
1349 tmumx(j) = tmumx(j) + mumx(j,k)
1350 tmumy(j) = tmumy(j) + mumy(j,k)
1351 tmumz(j) = tmumz(j) + mumz(j,k)
1352 END DO
1353 END DO
1354C
1355 RETURN

◆ wrrinf()

subroutine wrrinf ( character, dimension(*) title,
r,
integer n )

Definition at line 1471 of file hm_read_rbe3.F.

1472#include "implicit_f.inc"
1473#include "units_c.inc"
1474c !DECLARATIONS
1475 INTEGER N
1476 my_real
1477 . r(n)
1478 CHARACTER TITLE*(*)
1479C----------------------------
1480 INTEGER I
1481 write(iout, *)title,(r(i),i=1,n)
1482c print *,TITLE,(R(I),I=1,N)
1483 RETURN

◆ zero1()

subroutine zero1 ( r,
integer n )

Definition at line 1494 of file hm_read_rbe3.F.

1495#include "implicit_f.inc"
1496c !DECLARATIONS
1497 INTEGER N
1498 my_real
1499 . r(n)
1500C----------------------------
1501 INTEGER I
1502 DO i = 1,n
1503 r(i) = zero
1504 ENDDO
1505 RETURN