13 type(hecmwst_local_mesh),
intent(in) :: hecMESH
14 integer(kind=kint),
intent(in) :: numcomp
15 type(hecmwst_result_data),
intent(inout) :: fstrRESULT
19 call hecmw_nullify_result_data(fstrresult)
20 fstrresult%ng_component = 0
21 fstrresult%nn_component = numcomp
22 fstrresult%ne_component = 0
23 allocate( fstrresult%nn_dof(numcomp) )
24 allocate( fstrresult%node_label(numcomp) )
25 allocate( fstrresult%node_val_item(numcomp*hecmesh%n_dof*hecmesh%n_node))
30 type(hecmwst_result_data),
intent(inout) :: fstrRESULT
31 type(hecmwst_local_mesh),
intent(in) :: hecMESH
32 integer(kind=kint),
intent(in) :: comp_index
33 integer(kind=kint),
intent(in) :: ndof
34 character(len=HECMW_NAME_LEN),
intent(in) :: label
35 real(kind=kreal),
intent(in) :: vect(:)
37 integer(kind=kint) :: i, k, alldof, offset
40 fstrresult%nn_dof(comp_index) = ndof
41 fstrresult%node_label(comp_index) = label
42 alldof = fstrresult%nn_component*ndof
43 offset = ndof*(comp_index-1)
44 do i=1, hecmesh%n_node
46 fstrresult%node_val_item(alldof*(i-1) + k + offset) = vect(ndof*(i-1) + k)
69 fstrDYNAMIC, fstrRESULT, fstrPARAM, &
70 fstrCPL, fstrFREQ, hecLagMAT, restart_step_num)
74 type(hecmwst_local_mesh) :: hecMESH
75 type(hecmwst_matrix) :: hecMAT
78 type(hecmwst_result_data) :: fstrRESULT
83 type(hecmwst_matrix_lagrange) :: hecLagMAT
84 integer(kind=kint) :: restart_step_num
89 integer(kind=kint) :: numnode, numelm, startmode, endmode, nummode, ndof, im, in, ntotal, vistype
90 integer(kind=kint) :: numfreq, idnode, numdisp
91 integer(kind=kint) :: freqiout(3)
92 integer(kind=kint),
parameter :: ilogin = 9056
93 real(kind=kreal),
allocatable :: eigenvalue(:), loadvecre(:), loadvecim(:)
94 real(kind=kreal),
allocatable :: bjre(:), bjim(:), dvare(:), dvaim(:), disp(:), vel(:), acc(:)
95 real(kind=kreal),
allocatable :: dispre(:), dispim(:), velre(:), velim(:), accre(:), accim(:)
96 real(kind=kreal) :: freq, omega, val, dx, dy, dz, f_start, f_end
97 real(kind=kreal) :: t_start, t_end, time, dxi, dyi, dzi
100 numnode = hecmesh%nn_internal
101 numelm = hecmesh%n_elem
103 startmode = fstrfreq%start_mode
104 endmode = fstrfreq%end_mode
105 nummode = endmode - startmode +1
107 freqdata%numMode = nummode
108 freqdata%numNodeDOF = numnode*ndof
109 allocate(freqdata%eigOmega(nummode))
110 allocate(freqdata%eigVector(numnode*ndof, nummode))
112 call setupfreqparam(fstrdynamic, f_start, f_end, numfreq, freqdata%rayAlpha, freqdata%rayBeta, idnode, vistype, freqiout)
113 write(*,*)
"Rayleigh alpha:", freqdata%rayAlpha
114 write(*,*)
"Rayleigh beta:", freqdata%rayBeta
115 write(
ilog,*)
"Rayleigh alpha:", freqdata%rayAlpha
116 write(
ilog,*)
"Rayleigh beta:", freqdata%rayBeta
119 allocate(eigenvalue(nummode))
120 allocate(loadvecre(numnode*ndof))
121 allocate(loadvecim(numnode*ndof))
122 allocate(bjre(nummode))
123 allocate(bjim(nummode))
124 allocate(dvare(numnode*ndof))
125 allocate(dvaim(numnode*ndof))
126 allocate(disp(numnode*ndof))
127 allocate(vel(numnode*ndof))
128 allocate(acc(numnode*ndof))
129 allocate(dispre(numnode*ndof))
130 allocate(dispim(numnode*ndof))
131 allocate(velre(numnode*ndof))
132 allocate(velim(numnode*ndof))
133 allocate(accre(numnode*ndof))
134 allocate(accim(numnode*ndof))
139 write(*,*)
"--frequency analysis--"
140 write(*, *)
"read from=", trim(fstrfreq%eigenlog_filename)
141 write(
ilog,*)
"read from=", trim(fstrfreq%eigenlog_filename)
142 write(*, *)
"start mode=", startmode
143 write(
ilog,*)
"start mode=", startmode
144 write(*, *)
"end mode=", endmode
145 write(
ilog,*)
"end mode=", endmode
146 open(unit=ilogin, file=trim(fstrfreq%eigenlog_filename), status=
"OLD", action=
"READ")
147 call read_eigen_values(ilogin, startmode, endmode, eigenvalue, freqdata%eigOmega)
156 write(*,*)
"calc mass matrix"
157 call calcmassmatrix(fstrparam, hecmesh, hecmat, fstrsolid, fstreig, heclagmat)
158 write(*,*)
"scale eigenvector"
161 write(*, *)
"start frequency:", f_start
162 write(
ilog,*)
"start frequency:", f_start
163 write(*, *)
"end frequency:", f_end
164 write(
ilog,*)
"end frequency:", f_end
165 write(* ,*)
"number of the sampling points", numfreq
166 write(
ilog,*)
"number of the sampling points", numfreq
167 write(* ,*)
"monitor nodeid=", idnode
168 write(
ilog,*)
"monitor nodeid=", idnode
171 freq = (f_end-f_start)/dble(numfreq)*dble(im) + f_start
172 omega = 2.0d0 * 3.14159265358979d0 * freq
174 call calcfreqcoeff(freqdata, loadvecre, loadvecim, omega, bjre, bjim)
177 dx = sqrt(dvare(3*(idnode-1)+1)**2 + dvaim(3*(idnode-1)+1)**2)
178 dy = sqrt(dvare(3*(idnode-1)+2)**2 + dvaim(3*(idnode-1)+2)**2)
179 dz = sqrt(dvare(3*(idnode-1)+3)**2 + dvaim(3*(idnode-1)+3)**2)
180 val = sqrt(dx**2 + dy**2 + dz**2)
181 write(*, *) freq,
"[Hz] : ", val
182 write(
ilog, *) freq,
"[Hz] : ", val
183 disp(:) = abs(cmplx(dvare(:), dvaim(:)))
185 call calcvelvector(freqdata, omega, bjre, bjim, dvare, dvaim)
186 vel(:) = abs(cmplx(dvare(:), dvaim(:)))
188 call calcaccvector(freqdata, omega, bjre, bjim, dvare, dvaim)
189 acc(:) = abs(cmplx(dvare(:), dvaim(:)))
192 write(*, *) freq,
"[Hz] : ", im,
".res"
193 write(
ilog,*) freq,
"[Hz] : ", im,
".res"
196 if(
ivisual==1 .and. vistype==1)
then
197 write(*, *) freq,
"[Hz] : ", im,
".vis"
198 write(
ilog,*) freq,
"[Hz] : ", im,
".vis"
204 write(*, *)
"start time:", t_start
205 write(
ilog,*)
"start time:", t_start
206 write(*, *)
"end time:", t_end
207 write(
ilog,*)
"end time:", t_end
208 write(*, *)
"frequency:", freq
209 write(
ilog,*)
"frequency:", freq
210 write(*, *)
"node id:", idnode
211 write(
ilog,*)
"node id:", idnode
212 write(*, *)
"num disp:", numdisp
213 write(
ilog,*)
"num disp:", numdisp
215 omega = 2.0d0 * 3.14159265358979d0 * freq
216 call calcfreqcoeff(freqdata, loadvecre, loadvecim, omega, bjre, bjim)
220 time = (t_end-t_start)/dble(numdisp)*dble(im-1) + t_start
222 dx = dvare(3*(idnode-1)+1)
223 dy = dvare(3*(idnode-1)+2)
224 dz = dvare(3*(idnode-1)+3)
225 dxi = dvaim(3*(idnode-1)+1)
226 dyi = dvaim(3*(idnode-1)+2)
227 dzi = dvaim(3*(idnode-1)+3)
232 write(*, *)
"time=", time,
" : ", im,
".res"
233 write(
ilog,*)
"time=", time,
" : ", im,
".res"
234 call outputdyna_resfile(hecmesh, time, im, dvare, dvaim, velre, velim, accre, accim, freqiout)
236 if(
ivisual==1 .and. vistype==2)
then
237 write(*, *)
"time=", time,
" : ", im,
".vis"
238 write(
ilog,*)
"time=", time,
" : ", im,
".vis"
239 call outputdyna_visfile(hecmesh, im, dvare, dvaim, velre, velim, accre, accim, freqiout)
243 deallocate(freqdata%eigOmega)
244 deallocate(freqdata%eigVector)
245 deallocate(eigenvalue)
246 deallocate(loadvecre)
247 deallocate(loadvecim)
266 integer(kind=kint),
intent(in) :: logfile
267 integer(kind=kint),
intent(in) :: startmode
268 integer(kind=kint),
intent(in) :: endmode
269 real(kind=kreal),
intent(inout) :: eigenvalue(:)
270 real(kind=kreal),
intent(inout) :: anglfreq(:)
272 integer(kind=kint) :: im, endflag, id
273 character(len=HECMW_MSG_LEN) :: line
274 real(kind=kreal) :: freq
281 read(logfile,
'(A80)', err=119) line
282 if(trim(adjustl(line)) ==
"NO. EIGENVALUE FREQUENCY (HZ) X Y Z X")
then
287 read(logfile,
'(A80)') line
290 read(logfile,
'(A80)') line
292 do im=1, (endmode-startmode+1)
293 read(logfile,
'(i5,3e12.4,a)', err=119) id, eigenvalue(im), anglfreq(im), freq, line
298 119
write(*,*)
"Error to find eigenvalue information from logfile"
299 write(
ilog,*)
"Error to find eigenvalue information from logfile"
305 integer(kind=kint),
intent(in) :: logfile
306 integer(kind=kint),
intent(in) :: startmode
307 integer(kind=kint),
intent(in) :: endmode
308 integer(kind=kint),
intent(in) :: numdof
309 integer(kind=kint),
intent(in) :: numnode
310 real(kind=kreal),
intent(inout) :: eigenvector(:, :)
312 integer(kind=kint) :: im, in, gblid, j, idx
313 real(kind=kreal) :: vec(6)
314 character(len=HECMW_MSG_LEN) :: line
321 read(logfile,
'(a80)', err=119,
end=119) line
322 if(line(1:9) ==
" Mode No.")
then
329 do im=1, (endmode-startmode+1)
332 read(logfile,
'(a80)', err=119,
end=119) line
333 if(line(1:9) ==
" Mode No.")
then
338 read(logfile,
'(a80)', err=119) line
339 read(logfile,
'(a80)', err=119) line
340 read(logfile,
'(a80)', err=119) line
341 read(logfile,
'(a80)', err=119) line
346 read(logfile,
'(i10,2e12.4)', err=119) gblid, (vec(j), j=1,2)
350 read(logfile,
'(i10,3e16.8)', err=119) gblid, (vec(j), j=1,3)
353 read(logfile,
'(i10,6e12.4)', err=119) gblid, (vec(j), j=1,6)
360 idx = (in-1)*numdof + j
361 eigenvector(idx,im) = vec(j)
368 119
write(*,*)
"Error to find eigenvector from logfile"
369 write(ilog,*)
"Error to find eigenvector from logfile"
376 type(hecmwst_local_mesh),
intent(in) :: hecMESH
377 integer(kind=kint),
intent(in) :: startmode
378 integer(kind=kint),
intent(in) :: endmode
379 integer(kind=kint),
intent(in) :: numdof
380 integer(kind=kint),
intent(in) :: numnode
381 real(kind=kreal),
intent(inout) :: eigenvector(:, :)
383 integer(kind=kint),
parameter :: compidx = 1
384 integer(kind=kint) :: imode, idx, ind, a, b, nallcomp, j
385 type(hecmwst_result_data) :: eigenres
386 character(len=HECMW_NAME_LEN) :: name
390 do imode=startmode, endmode
392 call hecmw_result_read_by_name(hecmesh, name, imode, eigenres)
395 do ind=1,eigenres%nn_component
396 nallcomp = nallcomp + eigenres%nn_dof(ind)
399 idx = imode - startmode + 1
402 a = (ind-1)*nallcomp + j
403 b = (ind-1)*numdof + j
404 eigenvector(b,imode) = eigenres%node_val_item(a)
414 type(hecmwst_result_data) :: res
417 if(
associated(res%nn_dof))
deallocate(res%nn_dof)
418 if(
associated(res%ne_dof))
deallocate(res%ne_dof)
419 if(
associated(res%node_label))
deallocate(res%node_label)
420 if(
associated(res%elem_label))
deallocate(res%elem_label)
421 if(
associated(res%node_val_item))
deallocate(res%node_val_item)
422 if(
associated(res%elem_val_item))
deallocate(res%elem_val_item)
427 type(hecmwst_result_data) :: res
432 nullify(res%node_label)
433 nullify(res%elem_label)
434 nullify(res%node_val_item)
435 nullify(res%elem_val_item)
443 type(hecmwst_local_mesh),
intent(in) :: hecMESH
444 real(kind=kreal),
intent(in) :: freq
445 integer(kind=kint),
intent(in) :: ifreq
446 real(kind=kreal),
intent(in) :: disp(:)
447 real(kind=kreal),
intent(in) :: vel(:)
448 real(kind=kreal),
intent(in) :: acc(:)
449 integer(kind=kint),
intent(in) :: iout(3)
451 integer(kind=kint) :: im
452 character(len=HECMW_HEADER_LEN) :: header
453 character(len=HECMW_MSG_LEN) :: comment
454 character(len=HECMW_NAME_LEN) :: label, nameid
455 real(kind=kreal) :: freqval(1)
460 comment=
'frequency_result'
461 call hecmw_result_init(hecmesh, ifreq, header, comment)
465 call hecmw_result_add(hecmw_result_dtype_global, 1, label, freqval)
467 if(iout(1) == 1)
then
469 call hecmw_result_add(hecmw_result_dtype_node, 3, label, disp)
471 if(iout(2) == 1)
then
473 call hecmw_result_add(hecmw_result_dtype_node, 3, label, vel)
475 if(iout(3) == 1)
then
477 call hecmw_result_add(hecmw_result_dtype_node, 3, label, acc)
479 call hecmw_result_write_by_name(nameid)
480 call hecmw_result_finalize()
486 type(hecmwst_local_mesh),
intent(in) :: hecmesh
487 integer(kind=kint),
intent(in) :: ifreq
488 real(kind=kreal),
intent(in) :: disp(:)
489 real(kind=kreal),
intent(in) :: vel(:)
490 real(kind=kreal),
intent(in) :: acc(:)
491 integer(kind=kint),
intent(in) :: iout(3)
493 type(hecmwst_result_data) :: fstrRESULT
494 character(len=HECMW_NAME_LEN) :: label
495 integer(kind=kint) :: ncomp, i
499 if(iout(i) == 1)
then
504 call fstr_freq_result_init(hecmesh, ncomp, fstrresult)
506 if(iout(1) == 1)
then
507 label =
'displace_abs'
508 call fstr_freq_result_add(fstrresult, hecmesh, ncomp, 3, label, disp)
511 if(iout(2) == 1)
then
512 label =
'velocity_abs'
513 call fstr_freq_result_add(fstrresult, hecmesh, ncomp, 3, label, vel)
516 if(iout(3) == 1)
then
517 label =
'acceleration_abs'
518 call fstr_freq_result_add(fstrresult, hecmesh, ncomp, 3, label, acc)
522 call fstr2hecmw_mesh_conv(hecmesh)
523 call hecmw_visualize_init
524 call hecmw_visualize( hecmesh, fstrresult, ifreq )
525 call hecmw2fstr_mesh_conv(hecmesh)
526 call hecmw_result_free(fstrresult)
531 type(hecmwst_local_mesh),
intent(in) :: hecMESH
532 type(fstr_freqanalysis),
intent(in) :: freqData
533 integer(kind=kint),
intent(in) :: numdof
534 real(kind=kreal),
intent(inout) :: loadvecre(:)
535 real(kind=kreal),
intent(inout) :: loadvecim(:)
537 integer(kind=kint),
parameter :: MAXNODE = 100
538 integer(kind=kint) :: sgrpID, is, ie, ic, nsurf, ic_type, outtype, node_index(MAXNODE)
539 integer(kind=kint) :: nn, iss, nodeid, dof_index, ndof
540 integer(kind=kint) :: i, j, k, l, m, isn, nsize
541 integer(kind=kint) :: iwk(60), nodLOCAL(20)
542 real(kind=kreal) :: vect(60), xx(20), yy(20), zz(20), forcere(3), forceim(3)
546 do i=1,freqdata%FLOAD_ngrp_tot
547 if(freqdata%FLOAD_ngrp_TYPE(i) == kfloadtype_surf)
then
548 sgrpid = freqdata%FLOAD_ngrp_ID(i)
549 dof_index = freqdata%FLOAD_ngrp_DOF(i)
552 forcere(dof_index) = freqdata%FLOAD_ngrp_valre(i)
553 forceim(dof_index) = freqdata%FLOAD_ngrp_valim(i)
555 is = hecmesh%surf_group%grp_index(sgrpid-1) + 1
556 ie = hecmesh%surf_group%grp_index(sgrpid)
558 ic = hecmesh%surf_group%grp_item(2*j-1)
559 nsurf = hecmesh%surf_group%grp_item(2*j)
560 ic_type = hecmesh%elem_type(ic)
561 nn = hecmw_get_max_node(ic_type)
562 isn = hecmesh%elem_node_index(ic-1)
564 nodlocal(k) = hecmesh%elem_node_item(isn+k)
565 xx(k) = hecmesh%node(3*nodlocal(k)-2)
566 yy(k) = hecmesh%node(3*nodlocal(k)-1)
567 zz(k) = hecmesh%node(3*nodlocal(k) )
569 iwk(ndof*(k-1)+l) = ndof*(nodlocal(k)-1)+l
573 call dl_c3_freq(ic_type, nn, xx, yy, zz, nsurf, forcere, vect, nsize)
575 loadvecre(iwk(k)) = loadvecre(iwk(k)) + vect(k)
578 call dl_c3_freq(ic_type, nn, xx, yy, zz, nsurf, forceim, vect, nsize)
580 loadvecim(iwk(k)) = loadvecim(iwk(k)) + vect(k)
589 subroutine dl_c3_freq(ETYPE, NN, XX, YY, ZZ, LTYPE, force, VECT, nsize)
591 integer(kind=kint),
intent(in) :: ETYPE
592 integer(kind=kint),
intent(in) :: NN
593 integer(kind=kint),
intent(in) :: LTYPE
594 real(kind=kreal),
intent(in) :: xx(:)
595 real(kind=kreal),
intent(in) :: yy(:)
596 real(kind=kreal),
intent(in) :: zz(:)
597 real(kind=kreal),
intent(in) :: force(3)
598 real(kind=kreal),
intent(inout) :: vect(:)
599 integer(kind=kint),
intent(inout) :: nsize
601 integer(kind=kint),
parameter :: NDOF = 3
602 real(kind=kreal) :: wg
603 integer(kind=kint) :: NOD(NN)
604 real(kind=kreal) :: elecoord(3, nn), localcoord(3)
605 real(kind=kreal) :: h(nn)
606 integer(kind=kint) :: I, IG2, NSUR, SURTYPE
609 call getsubface( etype, ltype, surtype, nod )
610 nsur = getnumberofnodes( surtype )
613 elecoord(1,i)=xx(nod(i))
614 elecoord(2,i)=yy(nod(i))
615 elecoord(3,i)=zz(nod(i))
618 vect(1:nsize) = 0.0d0
619 do ig2=1,numofquadpoints( surtype )
620 call getquadpoint( surtype, ig2, localcoord(1:2) )
621 call getshapefunc( surtype, localcoord(1:2), h(1:nsur) )
623 wg=getweight( surtype, ig2 )
625 vect(3*nod(i)-2)=vect(3*nod(i)-2)+wg*h(i)*force(1)
626 vect(3*nod(i)-1)=vect(3*nod(i)-1)+wg*h(i)*force(2)
627 vect(3*nod(i) )=vect(3*nod(i) )+wg*h(i)*force(3)
634 type(hecmwst_local_mesh),
intent(in) :: hecMESH
635 type(fstr_freqanalysis),
intent(in) :: freqData
636 integer(kind=kint),
intent(in) :: numdof
637 real(kind=kreal),
intent(inout) :: loadvecre(:)
638 real(kind=kreal),
intent(inout) :: loadvecim(:)
640 integer(kind=kint) :: i, vecsize, ig, is, ie, in, nodeid, dof_index
644 do i=1, freqdata%FLOAD_ngrp_tot
645 if(freqdata%FLOAD_ngrp_TYPE(i) == kfloadtype_node)
then
646 ig = freqdata%FLOAD_ngrp_ID(i)
647 is = hecmesh%node_group%grp_index(ig-1) + 1
648 ie = hecmesh%node_group%grp_index(ig)
650 nodeid = hecmesh%node_group%grp_item(in)
651 dof_index = freqdata%FLOAD_ngrp_DOF(i)
652 loadvecre((nodeid-1)*numdof + dof_index) = loadvecre((nodeid-1)*numdof + dof_index) + freqdata%FLOAD_ngrp_valre(i)
653 loadvecim((nodeid-1)*numdof + dof_index) = loadvecim((nodeid-1)*numdof + dof_index) + freqdata%FLOAD_ngrp_valim(i)
660 subroutine calcmassmatrix(fstrPARAM, hecMESH, hecMAT, fstrSOLID, fstrEIG, hecLagMAT)
662 type(fstr_param),
intent(in) :: fstrPARAM
663 type(hecmwst_local_mesh),
intent(in) :: hecMESH
664 type(hecmwst_matrix),
intent(inout) :: hecMAT
665 type(fstr_solid),
intent(inout) :: fstrSOLID
666 type(fstr_eigen),
intent(inout) :: fstrEIG
667 type(hecmwst_matrix_lagrange),
intent(inout) :: hecLagMAT
669 integer(kind=kint) :: ntotal
673 fstrsolid%dunode = 0.d0
674 call fstr_stiffmatrix( hecmesh, hecmat, fstrsolid, 0.d0, 0.d0 )
675 call fstr_addbc(1, hecmesh, hecmat, fstrsolid, fstrparam, heclagmat, 2)
677 call setmass(fstrsolid, hecmesh, hecmat, fstreig)
683 type(fstr_eigen),
intent(in) :: fstrEIG
684 integer(kind=kint),
intent(in) :: ntotaldof
685 integer(kind=kint),
intent(in) :: nmode
686 real(kind=kreal),
intent(inout) :: eigenvector(:, :)
688 integer(kind=kint) :: imode, idof
689 real(kind=kreal) :: mas
695 mas = mas + fstreig%mass(idof)*eigenvector(idof,imode)**2
698 eigenvector(idof,imode) = eigenvector(idof,imode) / sqrt(mas)
705 type(fstr_eigen),
intent(in) :: fstreig
706 real(kind=kreal),
intent(in) :: eigenvector(:, :)
707 integer(kind=kint),
intent(in) :: imode
708 integer(kind=kint),
intent(in) :: jmode
709 real(kind=kreal),
intent(inout) :: prod
711 integer(kind=kint) :: idof, s
713 s =
size(eigenvector(:,1))
717 prod = prod + eigenvector(idof,imode)*fstreig%mass(idof)*eigenvector(idof,jmode)
724 integer(kind=kint),
intent(in) :: im
725 real(kind=kreal),
intent(in) :: vector(:)
727 integer(kind=kint) :: i, s
732 write(*,
'("eigenvec",i2.2,":[",e12.5", ")') im, vector(i)
734 write(*,
'(e12.5,", ")') vector(i)
736 write(*,
'(e12.5,"];")') vector(i)
745 real(kind=kreal),
intent(in) :: a(:)
746 real(kind=kreal),
intent(in) :: b(:)
747 real(kind=kreal),
intent(inout) :: c
750 c = dot_product(a, b)
757 type(fstr_freqanalysis_data),
intent(in) :: freqdata
758 real(kind=kreal),
intent(in) :: loadre(:)
759 real(kind=kreal),
intent(in) :: loadim(:)
760 real(kind=kreal),
intent(in) :: inpomega
761 real(kind=kreal),
intent(inout) :: bjre(:)
762 real(kind=kreal),
intent(inout) :: bjim(:)
764 integer(kind=kint) :: imode
765 real(kind=kreal) :: ujfr, ujfi, a, b, alp, beta
768 alp = freqdata%rayAlpha
769 beta = freqdata%rayBeta
771 do imode=1, freqdata%numMode
775 a = ujfr*(freqdata%eigOmega(imode)**2 - inpomega**2) + ujfi*(alp + beta*freqdata%eigOmega(imode)**2)*inpomega
776 b = (freqdata%eigOmega(imode)**2 - inpomega**2)**2 + ((alp + beta*freqdata%eigOmega(imode)**2)*inpomega)**2
779 a = ujfi*(freqdata%eigOmega(imode)**2 -inpomega**2) - ujfr*(alp + beta*freqdata%eigOmega(imode)**2)*inpomega
780 b = (freqdata%eigOmega(imode)**2 - inpomega**2)**2 + ((alp + beta*freqdata%eigOmega(imode)**2)*inpomega)**2
788 type(fstr_freqanalysis_data),
intent(in) :: freqdata
789 real(kind=kreal),
intent(in) :: bjre(:)
790 real(kind=kreal),
intent(in) :: bjim(:)
791 real(kind=kreal),
intent(inout) :: dispre(:)
792 real(kind=kreal),
intent(inout) :: dispim(:)
794 integer(kind=kint) :: imode
800 do imode=1, freqdata%numMode
801 dispre(:) = dispre(:) + bjre(imode)*freqdata%eigVector(:,imode)
802 dispim(:) = dispim(:) + bjim(imode)*freqdata%eigVector(:,imode)
809 type(fstr_freqanalysis_data),
intent(in) :: freqData
810 real(kind=kreal),
intent(in) :: omega
811 real(kind=kreal),
intent(in) :: bjre(:)
812 real(kind=kreal),
intent(in) :: bjim(:)
813 real(kind=kreal),
intent(inout) :: velre(:)
814 real(kind=kreal),
intent(inout) :: velim(:)
816 integer(kind=kint) :: imode
822 do imode=1, freqdata%numMode
823 velre(:) = velre(:) - omega * bjim(imode) * freqdata%eigVector(:,imode)
824 velim(:) = velim(:) + omega * bjre(imode) * freqdata%eigVector(:,imode)
830 type(fstr_freqanalysis_data),
intent(in) :: freqData
831 real(kind=kreal),
intent(in) :: omega
832 real(kind=kreal),
intent(in) :: bjre(:)
833 real(kind=kreal),
intent(in) :: bjim(:)
834 real(kind=kreal),
intent(inout) :: accre(:)
835 real(kind=kreal),
intent(inout) :: accim(:)
837 integer(kind=kint) :: imode
843 do imode=1, freqdata%numMode
844 accre(:) = accre(:) - omega**2 * bjre(imode) * freqdata%eigVector(:,imode)
845 accim(:) = accim(:) - omega**2 * bjim(imode) * freqdata%eigVector(:,imode)
850 subroutine setupfreqparam(fstrDYNAMIC, f_start, f_end, numfreq, raym, rayk, idnode, vistype, ioutl)
852 type(fstr_dynamic),
intent(in) :: fstrDYNAMIC
853 real(kind=kreal),
intent(inout) :: f_start
854 real(kind=kreal),
intent(inout) :: f_end
855 integer(kind=kint),
intent(inout) :: numfreq
856 real(kind=kreal),
intent(inout) :: raym
857 real(kind=kreal),
intent(inout) :: rayk
858 integer(kind=kint),
intent(inout) :: idnode
859 integer(kind=kint),
intent(inout) :: vistype
860 integer(kind=kint),
intent(inout) :: ioutl(3)
864 f_start = fstrdynamic%t_start
865 f_end = fstrdynamic%t_end
866 numfreq = fstrdynamic%n_step
867 raym = fstrdynamic%ray_m
868 rayk = fstrdynamic%ray_k
869 idnode = fstrdynamic%nout_monit
870 vistype = fstrdynamic%ngrp_monit
871 ioutl(1:3) = fstrdynamic%iout_list(1:3)
877 type(fstr_freqanalysis_data),
intent(in) :: freqData
878 real(kind=kreal),
intent(in) :: time
879 real(kind=kreal),
intent(in) :: omega
880 real(kind=kreal),
intent(in) :: bjre(:)
881 real(kind=kreal),
intent(in) :: bjim(:)
882 real(kind=kreal),
intent(inout) :: dispre(:)
883 real(kind=kreal),
intent(inout) :: dispim(:)
885 integer(kind=kint) :: imode, idf, s
886 complex(kind=kreal) :: a, b, c
891 a = exp(cmplx(0.0d0, omega*time))
893 do imode=1, freqdata%numMode
894 s =
size(freqdata%eigvector(:,imode))
895 b = cmplx(bjre(imode), bjim(imode)) * a
897 c = b*cmplx(freqdata%eigVector(idf,imode), 0.0d0)
898 dispre(idf) = dispre(idf) + dble(c)
899 dispim(idf) = dispim(idf) + imag(c)
907 type(fstr_freqanalysis_data),
intent(in) :: freqData
908 real(kind=kreal),
intent(in) :: time
909 real(kind=kreal),
intent(in) :: omega
910 real(kind=kreal),
intent(in) :: bjre(:)
911 real(kind=kreal),
intent(in) :: bjim(:)
912 real(kind=kreal),
intent(inout) :: velre(:)
913 real(kind=kreal),
intent(inout) :: velim(:)
915 integer(kind=kint) :: imode, idf, s
916 complex(kind=kreal) :: a, b, c
921 a = cmplx(0.0d0, 1.0d0)*cmplx(omega, 0.0d0)*exp(cmplx(0.0d0, omega*time))
923 do imode=1, freqdata%numMode
924 s =
size(freqdata%eigvector(:,imode))
925 b = cmplx(bjre(imode), bjim(imode)) * a
927 c = b*cmplx(freqdata%eigVector(idf,imode), 0.0d0)
928 velre(idf) = velre(idf) + dble(c)
929 velim(idf) = velim(idf) + imag(c)
937 type(fstr_freqanalysis_data),
intent(in) :: freqData
938 real(kind=kreal),
intent(in) :: time
939 real(kind=kreal),
intent(in) :: omega
940 real(kind=kreal),
intent(in) :: bjre(:)
941 real(kind=kreal),
intent(in) :: bjim(:)
942 real(kind=kreal),
intent(inout) :: accre(:)
943 real(kind=kreal),
intent(inout) :: accim(:)
945 integer(kind=kint) :: imode, idf, s
946 complex(kind=kreal) :: a, b, c
951 a = cmplx(-1.0d0, 0.0d0)*cmplx(omega**2, 0.0d0)*exp(cmplx(0.0d0, omega*time))
953 do imode=1, freqdata%numMode
954 s =
size(freqdata%eigvector(:,imode))
955 b = cmplx(bjre(imode), bjim(imode)) * a
957 c = b*cmplx(freqdata%eigVector(idf,imode), 0.0d0)
958 accre(idf) = accre(idf) + dble(c)
959 accim(idf) = accim(idf) + imag(c)
967 type(fstr_dynamic),
intent(in) :: fstrDYNAMIC
968 real(kind=kreal),
intent(inout) :: t_start
969 real(kind=kreal),
intent(inout) :: t_end
970 real(kind=kreal),
intent(inout) :: dynafreq
971 integer(kind=kint),
intent(inout) :: numdisp
974 t_start = fstrdynamic%gamma
975 t_end = fstrdynamic%beta
976 dynafreq = fstrdynamic%t_delta
977 numdisp = fstrdynamic%nout
983 type(hecmwst_local_mesh),
intent(in) :: hecMESH
984 real(kind=kreal),
intent(in) :: time
985 integer(kind=kint),
intent(in) :: istp
986 real(kind=kreal),
intent(in) :: dispre(:)
987 real(kind=kreal),
intent(in) :: dispim(:)
988 real(kind=kreal),
intent(in) :: velre(:)
989 real(kind=kreal),
intent(in) :: velim(:)
990 real(kind=kreal),
intent(in) :: accre(:)
991 real(kind=kreal),
intent(in) :: accim(:)
992 integer(kind=kint),
intent(in) :: iout(3)
994 integer(kind=kint) :: im, s
995 character(len=HECMW_HEADER_LEN) :: header
996 character(len=HECMW_MSG_LEN) :: comment
997 character(len=HECMW_NAME_LEN) :: label, nameid
998 real(kind=kreal),
allocatable :: absval(:)
1005 header=
'*fstrresult'
1006 comment=
'frequency_result'
1008 call hecmw_result_init(hecmesh, istp, header, comment)
1012 call hecmw_result_add(hecmw_result_dtype_global, 1, label, absval)
1014 if(iout(1) == 1)
then
1015 label=
'displacement_real'
1016 call hecmw_result_add(hecmw_result_dtype_node, 3, label, dispre)
1017 label=
'displacement_imag'
1018 call hecmw_result_add(hecmw_result_dtype_node, 3, label, dispim)
1019 label=
'displacement_abs'
1020 absval(:) = abs(cmplx(dispre(:), dispim(:)))
1021 call hecmw_result_add(hecmw_result_dtype_node, 3, label, absval)
1024 if(iout(2) == 1)
then
1025 label=
'velocity_real'
1026 call hecmw_result_add(hecmw_result_dtype_node, 3, label, velre)
1027 label=
'velocity_imag'
1028 call hecmw_result_add(hecmw_result_dtype_node, 3, label, velim)
1029 label=
'velocity_abs'
1030 absval(:) = abs(cmplx(velre(:), velim(:)))
1031 call hecmw_result_add(hecmw_result_dtype_node, 3, label, absval)
1034 if(iout(3) == 1)
then
1035 label=
'acceleration_real'
1036 call hecmw_result_add(hecmw_result_dtype_node, 3, label, accre)
1037 label=
'acceleration_imag'
1038 call hecmw_result_add(hecmw_result_dtype_node, 3, label, accim)
1039 label=
'acceleration_abs'
1040 absval(:) = abs(cmplx(velre(:), velim(:)))
1041 call hecmw_result_add(hecmw_result_dtype_node, 3, label, absval)
1044 call hecmw_result_write_by_name(nameid)
1045 call hecmw_result_finalize()
1053 type(hecmwst_local_mesh),
intent(inout) :: hecmesh
1054 integer(kind=kint),
intent(in) :: istp
1055 real(kind=kreal),
intent(in) :: dispre(:)
1056 real(kind=kreal),
intent(in) :: dispim(:)
1057 real(kind=kreal),
intent(in) :: velre(:)
1058 real(kind=kreal),
intent(in) :: velim(:)
1059 real(kind=kreal),
intent(in) :: accre(:)
1060 real(kind=kreal),
intent(in) :: accim(:)
1061 integer(kind=kint),
intent(in) :: iout(3)
1063 type(hecmwst_result_data) :: fstrRESULT
1064 character(len=HECMW_NAME_LEN) :: label
1065 integer(kind=kint) :: s, ncomp, i
1066 real(kind=kreal),
allocatable :: absval(:)
1074 if(iout(i) == 1)
then
1079 call fstr_freq_result_init(hecmesh, ncomp, fstrresult)
1083 if(iout(1) == 1)
then
1084 label =
'displace_real'
1085 call fstr_freq_result_add(fstrresult, hecmesh, ncomp, 3, label, dispre)
1088 label =
'displace_imag'
1089 call fstr_freq_result_add(fstrresult, hecmesh, ncomp, 3, label, dispim)
1092 label =
'displace_abs'
1093 absval(:) = abs(cmplx(dispre(:), dispim(:)))
1094 call fstr_freq_result_add(fstrresult, hecmesh, ncomp, 3, label, absval)
1098 if(iout(2) == 1)
then
1099 label =
'velocity_real'
1100 call fstr_freq_result_add(fstrresult, hecmesh, ncomp, 3, label, velre)
1103 label =
'velocity_imag'
1104 call fstr_freq_result_add(fstrresult, hecmesh, ncomp, 3, label, velim)
1107 label =
'velocity_abs'
1108 absval(:) = abs(cmplx(velre(:), velim(:)))
1109 call fstr_freq_result_add(fstrresult, hecmesh, ncomp, 3, label, absval)
1113 if(iout(3) == 1)
then
1114 label =
'acceleration_real'
1115 call fstr_freq_result_add(fstrresult, hecmesh, ncomp, 3, label, accre)
1118 label =
'acceleration_imag'
1119 call fstr_freq_result_add(fstrresult, hecmesh, ncomp, 3, label, accim)
1122 label =
'acceleration_abs'
1123 absval(:) = abs(cmplx(accre(:), accim(:)))
1124 call fstr_freq_result_add(fstrresult, hecmesh, ncomp, 3, label, absval)
1128 call fstr2hecmw_mesh_conv(hecmesh)
1129 call hecmw_visualize_init
1130 call hecmw_visualize( hecmesh, fstrresult, istp )
1131 call hecmw2fstr_mesh_conv(hecmesh)
1132 call hecmw_result_free(fstrresult)
subroutine nullify_result_data(res)
subroutine free_result_data(res)
subroutine calcdispvector(freqData, bjRe, bjIm, dispRe, dispIm)
subroutine read_eigen_vector(logfile, startmode, endmode, numdof, numnode, eigenvector)
subroutine calcdotproduct(a, b, c)
subroutine output_resfile(hecMESH, freq, ifreq, disp, vel, acc, iout)
subroutine fstr_solve_frequency_analysis(hecMESH, hecMAT, fstrSOLID, fstrEIG, fstrDYNAMIC, fstrRESULT, fstrPARAM, fstrCPL, fstrFREQ, hecLagMAT, restart_step_num)
subroutine read_eigen_vector_res(hecMESH, startmode, endmode, numdof, numnode, eigenvector)
subroutine calcaccvector(freqData, omega, bjRe, bjIm, accRe, accIm)
subroutine extract_surf2node(hecMESH, freqData, numdof, loadvecRe, loadvecIm)
subroutine writeoutvector(im, vector)
subroutine calcdispvectortime(freqData, time, omega, bjRe, bjIm, dispRe, dispIm)
subroutine checkorthvector(fstrEIG, eigenvector, imode, jmode, prod)
subroutine calcfreqcoeff(freqData, loadRe, loadIm, inpOmega, bjRe, bjIm)
subroutine scaleeigenvector(fstrEIG, ntotaldof, nmode, eigenvector)
subroutine output_visfile(hecMESH, ifreq, disp, vel, acc, iout)
subroutine setupdynaparam(fstrDYNAMIC, t_start, t_end, dynafreq, numdisp)
subroutine calcvelvectortime(freqData, time, omega, bjRe, bjIm, velRe, velIm)
subroutine setupfreqparam(fstrDYNAMIC, f_start, f_end, numfreq, raym, rayk, idnode, vistype, ioutl)
subroutine outputdyna_resfile(hecMESH, time, istp, dispre, dispim, velre, velim, accre, accim, iout)
subroutine read_eigen_values(logfile, startmode, endmode, eigenvalue, anglfreq)
subroutine calcmassmatrix(fstrPARAM, hecMESH, hecMAT, fstrSOLID, fstrEIG, hecLagMAT)
subroutine assemble_nodeload(hecMESH, freqData, numdof, loadvecRe, loadvecIm)
subroutine calcaccvectortime(freqData, time, omega, bjRe, bjIm, accRe, accIm)
subroutine outputdyna_visfile(hecMESH, istp, dispre, dispim, velre, velim, accre, accim, iout)
subroutine dl_c3_freq(ETYPE, NN, XX, YY, ZZ, LTYPE, force, VECT, nsize)
subroutine calcvelvector(freqData, omega, bjRe, bjIm, velRe, velIm)
This module contains steady state frequency analysis.
subroutine fstr_freq_result_init(hecMESH, numcomp, fstrRESULT)
subroutine fstr_freq_result_add(fstrRESULT, hecMESH, comp_index, ndof, label, vect)
This module provides a function to deal with prescribed displacement.
Set up lumped mass matrix.
This module provides function to calculate tangent stiffness matrix.
This module defines common data and basic structures for analysis.
integer(kind=kint), pointer iresult
integer(kind=kint), parameter ilog
FILE HANDLER.
integer(kind=kint), pointer ivisual
HECMW to FSTR Mesh Data Converter. Converting Connectivity of Element Type 232, 342 and 352.
Data for coupling analysis.
Data for DYNAMIC ANSLYSIS (fstrDYNAMIC)
Package of data used by Lanczos eigenvalue solver.
FSTR INNER CONTROL PARAMETERS (fstrPARAM)