OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
hm_read_damp.F File Reference
#include "implicit_f.inc"
#include "com01_c.inc"
#include "com04_c.inc"
#include "param_c.inc"
#include "units_c.inc"
#include "sphcom.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine hm_read_damp (dampr, igrnod, iskn, lsubmodel, unitab, snpc1, npc1, ndamp_vrel_rby, igrpart, damp_range_part)

Function/Subroutine Documentation

◆ hm_read_damp()

subroutine hm_read_damp ( dampr,
type (group_), dimension(ngrnod), target igrnod,
integer, dimension(liskn,*) iskn,
type(submodel_data), dimension(*), intent(in) lsubmodel,
type (unit_type_), intent(in) unitab,
integer, intent(in) snpc1,
integer, dimension(snpc1), intent(in) npc1,
integer, intent(inout) ndamp_vrel_rby,
type (group_), dimension(ngrpart) igrpart,
integer, dimension(npart), intent(inout) damp_range_part )
Parameters
[in,out]damp_range_partflag to compute the damping range

Definition at line 42 of file hm_read_damp.F.

44C-----------------------------------------------
45C M o d u l e s
46C-----------------------------------------------
47 USE message_mod
48 USE groupdef_mod
49 USE submodel_mod
51 USE unitab_mod
53 USE damping_range_compute_param_mod
54C-----------------------------------------------
55C I m p l i c i t T y p e s
56C-----------------------------------------------
57#include "implicit_f.inc"
58C-----------------------------------------------
59C C o m m o n B l o c k s
60C-----------------------------------------------
61#include "com01_c.inc"
62#include "com04_c.inc"
63#include "param_c.inc"
64#include "units_c.inc"
65#include "sphcom.inc"
66C-----------------------------------------------
67C D u m m y A r g u m e n t s
68C-----------------------------------------------
69 TYPE (UNIT_TYPE_),INTENT(IN) ::UNITAB
70 TYPE(SUBMODEL_DATA),INTENT(IN)::LSUBMODEL(*)
71 INTEGER ISKN(LISKN,*)
72 INTEGER, INTENT(IN) :: SNPC1,NPC1(SNPC1)
73 INTEGER, INTENT(INOUT) :: NDAMP_VREL_RBY
74 my_real dampr(nrdamp,*)
75 INTEGER, INTENT(INOUT) :: DAMP_RANGE_PART(NPART) !< flag to compute the damping range
76C-----------------------------------------------
77 TYPE (GROUP_) ,TARGET, DIMENSION(NGRNOD) :: IGRNOD
78 TYPE (GROUP_) ,DIMENSION(NGRPART) :: IGRPART
79C-----------------------------------------------
80C E x t e r n a l F u n c t i o n s
81C-----------------------------------------------
82 INTEGER NGR2USR
83 EXTERNAL ngr2usr
84C-----------------------------------------------
85C L o c a l V a r i a b l e s
86C-----------------------------------------------
87 INTEGER I,J,ID,JGRN,ISK,FL_VREL,FL_FREQ_RANGE,ITYPE
88 INTEGER NB_PAS,RANGE,FLINT,FLG_PRI,SUB_INDEX
89 INTEGER FUNC_ID,RBODY_ID,IFUN,IGR,GRPART,USR_GRPART
91 . factb,tstart,tstop,
92 . alpha,beta,alpha_y,beta_y,alpha_z,beta_z,alpha_xx,beta_xx,alpha_yy,
93 . beta_yy,alpha_zz,beta_zz,cdamp_mx,cdamp_my,cdamp_mz,
94 . dv2_mx,dv2_my,dv2_mz,freq,xscale,alpha_x,
95 . cdamp,freq_low,freq_high,maxwell_alpha(3),maxwell_tau(3)
96 CHARACTER(LEN=NCHARTITLE) :: TITR,KEY
97!
98 INTEGER, DIMENSION(:), POINTER :: INGR2USR
99 LOGICAL IS_AVAILABLE
100 LOGICAL FULL_FORMAT
101!-----------------------------------------------
102! DAMPR(1,I) : user_id
103! DAMPR(2,I) : IGRNOD ID
104! DAMPR(3,I) : alpha_x
105! DAMPR(4,I) : beta_x
106! DAMPR(5,I) : alpha_y
107! DAMPR(6,I) : beta_y
108! DAMPR(7,I) : alpha_z
109! DAMPR(8,I) : beta_z
110! DAMPR(9,I) : alpha_xx
111! DAMPR(10,I) : beta_xx
112! DAMPR(11,I) : alpha_yy
113! DAMPR(12,I) : beta_yy
114! DAMPR(13,I) : alpha_zz
115! DAMPR(14,I) : beta_zz
116! DAMPR(15,I) : ISK
117! DAMPR(16,I) : alpha of FUNC, FACTB now is fixed to 1
118! DAMPR(17,I) : TSTART ! to be initialized
119! DAMPR(18,I) : TSTOP ! to be initialized
120! DAMPR(19,I) : NB_PAS (0 if not INTER, 0 : using velocity in global system) ?
121! DAMPR(20,I) : RANGE (0 if not INTER)
122! DAMPR(21,I) : (0 INTER; 1 VREL; 0 FREQ) 2 for FUNC ?
123! Itype now -> 0:/DAMP; 1:/DAMP/INTER 2:/DAMP/VREL 3:/DAMP/FREQ 4:/DAMP/FUNCT
124! DAMPR(22,I) = Alpha2_x (only for VREL )
125! DAMPR(23,I) = Alpha2_y (only for VREL )
126! DAMPR(24,I) = Alpha2_z (only for VREL )
127! DAMPR(25,I) = RBODY_ID (only for VREL )
128! DAMPR(26,I) = IFUN (for VREL&FUNCT )
129! DAMPR(27,I) = XSCALE (only for VREL )
130! DAMPR(28,I) = FREQ (only for FREQ )
131! DAMPR(29,I) = ZERO (only for FREQ )
132! DAMPR(30,I) = ZERO (only for FREQ )
133! DAMPR(31,I) = (0 INTER; 0 VREL; 1 FREQ)
134! DAMPR(32:34,I) = MAXWELL_ALPHA(1:3) (only for FREQ ) alpha_x,alpha_y,alpha_z of FUNCT
135! DAMPR(35:37,I) = MAXWELL_TAU(1:3) (only for FREQ ) alpha_xx,alpha_yy,alpha_zz of FUNCT
136!======================================================================|
137 is_available = .false.
138 WRITE(iout,1000)
139C--------------------------------------------------
140C START BROWSING MODEL /DAMP
141C--------------------------------------------------
142 CALL hm_option_start('/DAMP')
143C--------------------------------------------------
144C BROWSING MODEL DAMP 1->NDAMP
145C--------------------------------------------------
146 DO i=1,ndamp
147C--------------------------------------------------
148C EXTRACT DATAS OF /DAMP/... LINE
149C--------------------------------------------------
150 CALL hm_option_read_key(lsubmodel,
151 . option_id = id,
152 . option_titr = titr,
153 . submodel_index = sub_index,
154 . keyword2=key)
155 full_format = .false.
156C--------------------------------------------------
157C HIDDEN FLAG FACTB
158C--------------------------------------------------
159C IF(NBLINES == 2) THEN
160C IREC=IREC+1
161C READ(IIN,REC=IREC,FMT=FMT_F) FACTB
162C ENDIF
163C--> SET TO 1.0
164C--------------------------------------------------
165 flint = 0
166 fl_vrel = 0
167 fl_freq_range = 0
168 itype = 0
169 factb = one
170C
171 IF(key(1:5)=='INTER')THEN
172 flint = 1
173 itype = 1
174 CALL hm_get_intv('Nb_time_step',nb_pas,is_available,lsubmodel)
175 CALL hm_get_intv('Range',range,is_available,lsubmodel)
176 CALL hm_get_intv('grnod_id',jgrn,is_available,lsubmodel)
177 CALL hm_get_intv('skew_id',isk,is_available,lsubmodel)
178 IF(isk == 0 .AND. sub_index /= 0 ) isk = lsubmodel(sub_index)%SKEW
179 CALL hm_get_boolv('Mass_Damp_Factor_Option',full_format,is_available)
180C--------------------------------------------------
181C EXTRACT DATAS (REAL VALUES)
182C--------------------------------------------------
183 CALL hm_get_floatv('Alpha',alpha,is_available,lsubmodel,unitab)
184 CALL hm_get_floatv('Beta',beta,is_available,lsubmodel,unitab)
185 CALL hm_get_floatv('Tstart',tstart,is_available,lsubmodel,unitab)
186 CALL hm_get_floatv('Tstop',tstop,is_available,lsubmodel,unitab)
187 CALL hm_get_floatv('Alpha_yy',alpha_yy,is_available,lsubmodel,unitab)
188 CALL hm_get_floatv('Beta_yy',beta_yy,is_available,lsubmodel,unitab)
189 CALL hm_get_floatv('Alpha_zz',alpha_zz,is_available,lsubmodel,unitab)
190 CALL hm_get_floatv('Beta_zz',beta_zz,is_available,lsubmodel,unitab)
191C--------------------------------------------------
192 IF (nb_pas == 0) nb_pas = 20
193 WRITE(iout,1300)
194 WRITE(iout,1400) nb_pas
195 WRITE(iout,1600) range
196 idamp_rdof = idamp_rdof+1
197 kcontact = 1
198 dampr(19,i) = nb_pas
199 dampr(20,i) = range
200 dampr(21,i) = 0
201 ELSEIF(key(1:4).EQ.'VREL')THEN
202 fl_vrel = 1
203 itype = 2
204C--------------------------------------------------
205C EXTRACT DATAS (INTEGER VALUES)
206C--------------------------------------------------
207 CALL hm_get_intv('grnod_id',jgrn,is_available,lsubmodel)
208 CALL hm_get_intv('skew_id',isk,is_available,lsubmodel)
209 IF(isk == 0 .AND. sub_index .NE. 0 ) isk = lsubmodel(sub_index)%SKEW
210 CALL hm_get_intv('RbodyID',rbody_id,is_available,lsubmodel)
211 CALL hm_get_intv('FuncID',func_id,is_available,lsubmodel)
212C--------------------------------------------------
213C EXTRACT DATAS (REAL VALUES)
214C--------------------------------------------------
215 CALL hm_get_floatv('Tstart',tstart,is_available,lsubmodel,unitab)
216 CALL hm_get_floatv('Tstop',tstop,is_available,lsubmodel,unitab)
217 CALL hm_get_floatv('Freq',freq,is_available,lsubmodel,unitab)
218 CALL hm_get_floatv('Xscale',xscale,is_available,lsubmodel,unitab)
219 CALL hm_get_floatv('Alpha_x',cdamp_mx,is_available,lsubmodel,unitab)
220 CALL hm_get_floatv('Alpha_y',cdamp_my,is_available,lsubmodel,unitab)
221 CALL hm_get_floatv('Alpha_z',cdamp_mz,is_available,lsubmodel,unitab)
222 CALL hm_get_floatv('Alpha2_x',dv2_mx,is_available,lsubmodel,unitab)
223 CALL hm_get_floatv('Alpha2_y',dv2_my,is_available,lsubmodel,unitab)
224 CALL hm_get_floatv('Alpha2_z',dv2_mz,is_available,lsubmodel,unitab)
225C--------------------------------------------------
226 IF (cdamp_my == zero) cdamp_my = cdamp_mx
227 IF (cdamp_mz == zero) cdamp_mz = cdamp_mx
228 IF (dv2_my == zero) dv2_my = dv2_mx
229 IF (dv2_mz == zero) dv2_mz = dv2_mx
230 IF (xscale == zero) xscale = one
231 alpha = zero
232 beta = zero
233C---------RBODY ----------------------------
234 IF (rbody_id /= 0) THEN
235 ndamp_vrel_rby = ndamp_vrel_rby + 1
236C RBODY merge - user ID of RBY is stored in DAMPR-> change to user id done after
237 ENDIF
238C---------FUNC ID user-----------------------------
239 ifun=0
240 IF (func_id /= 0) THEN
241 DO j=1,nfunct
242 IF (func_id == npc1(j)) THEN
243 ifun=j
244 EXIT
245 ENDIF
246 ENDDO
247 IF (ifun == 0)THEN ! Function not found
248 CALL ancmsg(msgid=3049,
249 . msgtype=msgerror,
250 . anmode=aninfo,
251 . i1=id,
252 . c1=titr,
253 . i2=func_id)
254 ENDIF
255 ENDIF
256C--------------------------------------------------
257 WRITE(iout,1700)
258 dampr(19,i) = 0
259 dampr(20,i) = 0
260 dampr(21,i) = 1
261 dampr(22,i) = dv2_mx
262 dampr(23,i) = dv2_my
263 dampr(24,i) = dv2_mz
264 dampr(25,i) = rbody_id
265 dampr(26,i) = ifun
266 dampr(27,i) = xscale
267 full_format = .true.
268C--------------------------------------------------
269 ELSEIF(key(1:4).EQ.'FREQ')THEN
270C--------------------------------------------------
271C Dapming in frequency range
272C--------------------------------------------------
273 itype = 3
274 fl_freq_range = 1
275C--------------------------------------------------
276C EXTRACT DATAS (INTEGER VALUES)
277C--------------------------------------------------
278 CALL hm_get_intv('grpart_id',grpart,is_available,lsubmodel)
279C---------Check part Id-----------------------------
280 IF(grpart/=0)THEN
281 igr = 0
282 DO j=1,ngrpart
283 IF (igrpart(j)%ID == grpart) THEN
284 igr=j
285 EXIT
286 END IF
287 END DO
288 IF(igr == 0) THEN
289 CALL ancmsg(msgid=3086,
290 . msgtype=msgerror,
291 . anmode=aninfo_blind_1,
292 . i1=id,
293 . c1=titr,
294 . i2=grpart)
295 ENDIF
296C---------Tag of the parts-----------------------------
297 DO j=1,igrpart(igr)%NENTITY
298 damp_range_part(igrpart(igr)%ENTITY(j)) = i
299 ENDDO
300 ELSE
301C---------Tag of all parts-----------------------------
302 DO j=1,npart
303 damp_range_part(j) = i
304 ENDDO
305 ENDIF
306 WRITE(iout,1900)
307 isk = 0
308 full_format = .true.
309 alpha = zero
310 beta = zero
311C--------------------------------------------------
312C EXTRACT DATAS (REAL VALUES)
313C--------------------------------------------------
314 CALL hm_get_floatv('Cdamp',cdamp,is_available,lsubmodel,unitab)
315 CALL hm_get_floatv('Tstart',tstart,is_available,lsubmodel,unitab)
316 CALL hm_get_floatv('Tstop',tstop,is_available,lsubmodel,unitab)
317 CALL hm_get_floatv('Freq_low',freq_low,is_available,lsubmodel,unitab)
318 CALL hm_get_floatv('Freq_high',freq_high,is_available,lsubmodel,unitab)
319C--------------------------------------------------
320 ELSEIF(key(1:5).EQ.'FUNCT')THEN
321C--------------------------------------------------
322 itype = 4
323C--------------------------------------------------
324C EXTRACT DATAS (INTEGER VALUES)
325C--------------------------------------------------
326 CALL hm_get_intv('grnod_id',jgrn,is_available,lsubmodel)
327 CALL hm_get_intv('FuncID',func_id,is_available,lsubmodel)
328C--------------------------------------------------
329C EXTRACT DATAS (REAL VALUES)
330C--------------------------------------------------
331 CALL hm_get_floatv('Alpha',alpha,is_available,lsubmodel,unitab)
332 CALL hm_get_floatv('Alpha_x',alpha_x,is_available,lsubmodel,unitab)
333 CALL hm_get_floatv('Alpha_y',alpha_y,is_available,lsubmodel,unitab)
334 CALL hm_get_floatv('Alpha_z',alpha_z,is_available,lsubmodel,unitab)
335 CALL hm_get_floatv('Alpha_xx',alpha_xx,is_available,lsubmodel,unitab)
336 CALL hm_get_floatv('Alpha_yy',alpha_yy,is_available,lsubmodel,unitab)
337 CALL hm_get_floatv('Alpha_zz',alpha_zz,is_available,lsubmodel,unitab)
338 isk = 0
339 tstart = zero
340 tstop=ep30
341 full_format = .true.
342 IF (alpha==zero) alpha = one ! default =1 to add in manual
343 factb = alpha
344 alpha = alpha_x
345!---------FUNC ID user-----------------------------
346 ifun=0
347 IF (func_id /= 0) THEN
348 DO j=1,nfunct
349 IF (func_id == npc1(j)) THEN
350 ifun=j
351 EXIT
352 ENDIF
353 ENDDO
354 IF (ifun == 0)THEN ! Function not found
355 CALL ancmsg(msgid=3049,
356 . msgtype=msgerror,
357 . anmode=aninfo,
358 . i1=id,
359 . c1=titr,
360 . i2=func_id)
361 ENDIF
362 ENDIF
363 dampr(4:nrdamp,i) = zero
364 WRITE(iout,2100)
365 ELSE
366C--------------------------------------------------
367C EXTRACT DATAS (INTEGER VALUES)
368C--------------------------------------------------
369 CALL hm_get_intv('grnod_id',jgrn,is_available,lsubmodel)
370 CALL hm_get_intv('skew_id',isk,is_available,lsubmodel)
371 IF(isk == 0 .AND. sub_index /= 0 ) isk = lsubmodel(sub_index)%SKEW
372 CALL hm_get_boolv('Mass_Damp_Factor_Option',full_format,is_available)
373C--------------------------------------------------
374C EXTRACT DATAS (REAL VALUES)
375C--------------------------------------------------
376 CALL hm_get_floatv('Alpha',alpha,is_available,lsubmodel,unitab)
377 CALL hm_get_floatv('Beta',beta,is_available,lsubmodel,unitab)
378 CALL hm_get_floatv('Tstart',tstart,is_available,lsubmodel,unitab)
379 CALL hm_get_floatv('Tstop',tstop,is_available,lsubmodel,unitab)
380 CALL hm_get_floatv('Alpha_y',alpha_y,is_available,lsubmodel,unitab)
381 CALL hm_get_floatv('Beta_y',beta_y,is_available,lsubmodel,unitab)
382 CALL hm_get_floatv('Alpha_z',alpha_z,is_available,lsubmodel,unitab)
383 CALL hm_get_floatv('Beta_z',beta_z,is_available,lsubmodel,unitab)
384 CALL hm_get_floatv('Alpha_xx',alpha_xx,is_available,lsubmodel,unitab)
385 CALL hm_get_floatv('Beta_xx',beta_xx,is_available,lsubmodel,unitab)
386 CALL hm_get_floatv('Alpha_yy',alpha_yy,is_available,lsubmodel,unitab)
387 CALL hm_get_floatv('Beta_yy',beta_yy,is_available,lsubmodel,unitab)
388 CALL hm_get_floatv('Alpha_zz',alpha_zz,is_available,lsubmodel,unitab)
389 CALL hm_get_floatv('Beta_zz',beta_zz,is_available,lsubmodel,unitab)
390C--------------------------------------------------
391 dampr(19,i) = 0
392 dampr(20,i) = 0
393 dampr(21,i) = 0
394C--------------------------------------------------
395 ENDIF ! IF(KEY(1:5)=='INTER')THEN
396C
397 DO j=0,numskw+min(1,nspcond)*numsph+nsubmod
398 IF(isk == iskn(4,j+1)) THEN
399 isk=j+1
400 GO TO 100
401 ENDIF
402 ENDDO
403 CALL ancmsg(msgid=137,anmode=aninfo,msgtype=msgerror,
404 . c1='DAMP',
405 . c2='DAMP',
406 . i1=id,i2=isk,c3=titr)
407 100 CONTINUE
408
409C
410 IF (tstop == zero) tstop=ep30
411C
412 dampr(1,i) = id
413 IF (fl_freq_range == 0) THEN
414 ingr2usr => igrnod(1:ngrnod)%ID
415 igr = ngr2usr(jgrn,ingr2usr,ngrnod)
416 IF (igr == 0) THEN
417 CALL ancmsg(msgid=171,
418 . msgtype=msgerror,
419 . anmode=aninfo,
420 . c1='RAYLEIGH DAMPING',
421 . i1= id,
422 . c2= titr,
423 . c3='NODE',
424 . i2=jgrn)
425 ENDIF
426 ENDIF
427 dampr(2,i) = igr
428 dampr(3,i) = alpha
429 dampr(4,i) = beta
430 dampr(15,i) = isk
431 dampr(17,i) = tstart
432 dampr(18,i) = tstop
433!
434 dampr(21,i) = itype
435C
436 IF (.NOT. full_format) THEN
437C-- reduced format only for /DAMP and /DAMP/INTER
438 dampr(5,i) = alpha
439 dampr(6,i) = beta
440 dampr(7,i) = alpha
441 dampr(8,i) = beta
442 dampr(9,i) = alpha
443 dampr(10,i) = beta
444 dampr(11,i) = alpha
445 dampr(12,i) = beta
446 dampr(13,i) = alpha
447 dampr(14,i) = beta
448 IF (flint==1) THEN
449 dampr(3,i) = zero
450 dampr(4,i) = zero
451 dampr(5,i) = zero
452 dampr(6,i) = zero
453 dampr(7,i) = zero
454 dampr(8,i) = zero
455 ENDIF
456 WRITE (iout,1100) jgrn,alpha,beta,factb,tstart,tstop
457 ELSE
458 SELECT CASE (itype)
459 CASE(0) !/DAMP
460 flg_pri = 1
461 dampr(3,i) = alpha
462 dampr(4,i) = beta
463 dampr(5,i) = alpha_y
464 dampr(6,i) = beta_y
465 dampr(7,i) = alpha_z
466 dampr(8,i) = beta_z
467 dampr(9,i) = alpha_xx
468 dampr(10,i) = beta_xx
469 dampr(11,i) = alpha_yy
470 dampr(12,i) = beta_yy
471 dampr(13,i) = alpha_zz
472 dampr(14,i) = beta_zz
473 WRITE (iout,1200) jgrn,iskn(4,isk),
474 . alpha,beta,alpha_y,beta_y,alpha_z,beta_z,
475 . alpha_xx,beta_xx,alpha_yy,beta_yy,alpha_zz,beta_zz,
476 . tstart,tstop
477 CASE(1) !/DAMP/INTER
478 dampr(3,i) = zero
479 dampr(4,i) = zero
480 dampr(5,i) = zero
481 dampr(6,i) = zero
482 dampr(7,i) = zero
483 dampr(8,i) = zero
484 dampr(9,i) = alpha
485 dampr(10,i) = beta
486 dampr(11,i) = alpha_yy
487 dampr(12,i) = beta_yy
488 dampr(13,i) = alpha_zz
489 dampr(14,i) = beta_zz
490 WRITE (iout,1500) jgrn,iskn(4,isk),
491 . alpha,beta,alpha_yy,beta_yy,
492 . alpha_zz,beta_zz,tstart,tstop
493 CASE(2) !/DAMP/VREL
494 dampr(3,i) = cdamp_mx
495 dampr(4,i) = zero
496 dampr(5,i) = cdamp_my
497 dampr(6,i) = zero
498 dampr(7,i) = cdamp_mz
499 dampr(8,i) = zero
500 dampr(9,i) = zero
501 dampr(10,i) = zero
502 dampr(11,i) = zero
503 dampr(12,i) = zero
504 dampr(13,i) = zero
505 dampr(14,i) = zero
506 WRITE (iout,1800) jgrn,iskn(4,isk),rbody_id,func_id,
507 . cdamp_mx,cdamp_my,cdamp_mz,
508 . dv2_mx,dv2_my,dv2_mz,
509 . freq,tstart,tstop
510 dampr(28,i) = freq
511 dampr(29,i) = zero
512 dampr(30,i) = zero
513 CASE(3) !/DAMP/FREQUENCY_RANGE
514 WRITE (iout,2000) grpart,cdamp,freq_low,freq_high,tstart,tstop
515C Automatic computation of parameters of the 3 maxwell components
516 CALL damping_range_compute_param(cdamp,freq_low,freq_high,maxwell_alpha,maxwell_tau)
517C
518 dampr(31,i) = one
519 dampr(32:34,i) = maxwell_alpha(1:3)
520 dampr(35:37,i) = maxwell_tau(1:3)
521 CASE(4) !/DAMP/FUNCT
522 alpha_x = alpha
523 dampr(3,i) = alpha_x
524 alpha = factb
525 dampr(5,i) = alpha_y
526 dampr(7,i) = alpha_z
527 dampr(9,i) = alpha_xx
528 dampr(11,i) = alpha_yy
529 dampr(13,i) = alpha_zz
530 dampr(26,i) = ifun ! take care of IFUN in split
531 dampr(32,i) = alpha_x
532 dampr(33,i) = alpha_y
533 dampr(34,i) = alpha_z
534 dampr(35,i) = alpha_xx
535 dampr(36,i) = alpha_yy
536 dampr(37,i) = alpha_zz
537 WRITE (iout,2200) jgrn,ifun,alpha,
538 . alpha_x,alpha_y,alpha_z,
539 . alpha_xx,alpha_yy,alpha_zz
540 END SELECT
541 END IF !(.NOT. FULL_FORMAT) THEN
542 dampr(16,i) = factb
543 END DO ! NDAMP
544C---
545 RETURN
546
547 1000 FORMAT(//
548 .' RAYLEIGH DAMPING '/
549 . ' ---------------------- ')
550 1100 FORMAT( 8x,'NODE GROUP ID . . . . . . . . .',i10
551 . /10x,'ALPHA. . . . . . . . . . . . . .',1pg20.13
552 . /10x,'BETA . . . . . . . . . . . . . .',1pg20.13
553 . /10x,'MAX TIME STEP FACTOR . . . . . .',1pg20.13
554 . /10x,'START TIME . . . . . . . . . . .',1pg20.13
555 . /10x,'STOP TIME . . . . . . . . . . .',1pg20.13)
556 1200 FORMAT( 10x,'NODE GROUP ID . . . . . . . . .',i10
557 . /10x,'SKEW ID . . . . . . . . . . .',i10
558 . /10x,'ALPHA IN X-DIRECTION. . . . . .',1pg20.13
559 . /10x,'BETA IN X-DIRECTION. . . . . .',1pg20.13
560 . /10x,'ALPHA IN Y-DIRECTION. . . . . .',1pg20.13
561 . /10x,'BETA IN Y-DIRECTION. . . . . .',1pg20.13
562 . /10x,'ALPHA IN Z-DIRECTION. . . . . .',1pg20.13
563 . /10x,'BETA IN Z-DIRECTION. . . . . .',1pg20.13
564 . /10x,'ALPHA IN RX-DIRECTION . . . . .',1pg20.13
565 . /10x,'BETA IN RX-DIRECTION . . . . .',1pg20.13
566 . /10x,'ALPHA IN RY-DIRECTION . . . . .',1pg20.13
567 . /10x,'BETA IN RY-DIRECTION . . . . .',1pg20.13
568 . /10x,'ALPHA IN RZ-DIRECTION . . . . .',1pg20.13
569 . /10x,'BETA IN RZ-DIRECTION . . . . .',1pg20.13
570 . /10x,'START TIME . . . . . . . . . . .',1pg20.13
571 . /10x,'STOP TIME . . . . . . . . . . .',1pg20.13)
572 1300 FORMAT(/,10x,'SELECTIVE RAYLEIGH DAMPING ON CONTACT NODES')
573 1400 FORMAT( 10x,'NUMBER OF TIME STEP . . . . . .',i10,/)
574 1500 FORMAT( 10x,'NODE GROUP ID . . . . . . . . .',i10
575 . /10x,'SKEW ID . . . . . . . . . . .',i10
576 . /10x,'ALPHA IN RX-DIRECTION . . . . .',1pg20.13
577 . /10x,'BETA IN RX-DIRECTION . . . . .',1pg20.13
578 . /10x,'ALPHA IN RY-DIRECTION . . . . .',1pg20.13
579 . /10x,'BETA IN RY-DIRECTION . . . . .',1pg20.13
580 . /10x,'ALPHA IN RZ-DIRECTION . . . . .',1pg20.13
581 . /10x,'BETA IN RZ-DIRECTION . . . . .',1pg20.13
582 . /10x,'START TIME . . . . . . . . . . .',1pg20.13
583 . /10x,'STOP TIME . . . . . . . . . . .',1pg20.13)
584 1600 FORMAT( 10x,'EXTENSION OF NODES SELECTION . ',i10,/)
585 1700 FORMAT(/,10x,'RAYLEIGH DAMPING WITH RELATIVE VELOCITIES')
586 1800 FORMAT( 10x,'NODE GROUP ID . . . . . . . . .',i10
587 . /10x,'SKEW ID . . . . . . . . . . . .',i10
588 . /10x,'RBODY ID . . . . . . . . . . . ',i10
589 . /10x,'DAMPING FUNCTION ID . . . . . .',i10
590 . /10x,'MASS DAMPING COEFFICIENT IN X-DIRECTION. . . . . .',1pg20.13
591 . /10x,'MASS DAMPING COEFFICIENT IN Y-DIRECTION. . . . . .',1pg20.13
592 . /10x,'MASS DAMPING COEFFICIENT IN Z-DIRECTION. . . . . .',1pg20.13
593 . /10x,'QUADRATIC MASS DAMPING COEFFICIENT IN X-DIRECTION.',1pg20.13
594 . /10x,'QUADRATIC MASS DAMPING COEFFICIENT IN Y-DIRECTION.',1pg20.13
595 . /10x,'QUADRATIC MASS DAMPING COEFFICIENT IN Z-DIRECTION.',1pg20.13
596 . /10x,'DAMPING FREQUENCY . . . . . . . . . . . . . . . . ',1pg20.13
597 . /10x,'START TIME . . . . . . . . . . . . . . . . . . . .',1pg20.13
598 . /10x,'STOP TIME . . . . . . . . . . . . . . . . . . . .',1pg20.13)
599 1900 FORMAT(/,10x,'DAMPING OVER FREQUENCY RANGE')
600 2000 FORMAT( 10x,'PART GROUP ID . . . . . . . . .',i10
601 . /10x,'DAMPING RATIO . . . . . . . . . . . . . . . . . . ',1pg20.13
602 . /10x,'LOWEST FREQUENCY . . . . . . . . . . . . . . . . .',1pg20.13
603 . /10x,'HIGHEST FREQUENCY. . . . . . . . . . . . . . . . .',1pg20.13
604 . /10x,'START TIME . . . . . . . . . . . . . . . . . . . .',1pg20.13
605 . /10x,'STOP TIME . . . . . . . . . . . . . . . . . . . .',1pg20.13)
606 2100 FORMAT(/,10x,'MASS DAMPING WITH INPUT FUNCTION')
607 2200 FORMAT( 10x,'NODE GROUP ID . . . . . . . . . . . . . . . . . .',i10
608 . /10x,'ALPHA FUNCTION ID . . . . . . . . . . . . . . . .',i10
609 . /10x,'ALPHA FUNCTION ORDINATE SCALE FACTOR . . . . . . ',1pg20.13
610 . /10x,'MASS DAMPING COEFFICIENT IN X-DIRECTION. . . . . ',1pg20.13
611 . /10x,'MASS DAMPING COEFFICIENT IN Y-DIRECTION. . . . . ',1pg20.13
612 . /10x,'MASS DAMPING COEFFICIENT IN Z-DIRECTION. . . . . ',1pg20.13
613 . /10x,'MASS DAMPING COEFFICIENT IN RX-DIRECTION. . . . .',1pg20.13
614 . /10x,'MASS DAMPING COEFFICIENT IN RY-DIRECTION. . . . .',1pg20.13
615 . /10x,'MASS DAMPING COEFFICIENT IN RZ-DIRECTION. . . . .',1pg20.13)
616C---
617 RETURN
#define my_real
Definition cppsort.cpp:32
#define alpha
Definition eval.h:35
subroutine hm_get_boolv(name, bval, is_available)
subroutine hm_get_floatv(name, rval, is_available, lsubmodel, unitab)
subroutine hm_get_intv(name, ival, is_available, lsubmodel)
subroutine hm_option_start(entity_type)
#define min(a, b)
Definition macros.h:20
initmumps id
integer, parameter nchartitle
integer nsubmod
integer function ngr2usr(iu, igr, ngr)
Definition nintrr.F:325
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:889