FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
m_fstr.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 ! Copyright (c) 2019 FrontISTR Commons
3 ! This software is released under the MIT License, see LICENSE.txt
4 !-------------------------------------------------------------------------------
5 ! If new header is supported, change according to following method. !
6 ! 1) Increase FSTR_CTRL_HEADER_NUMBER !
7 ! 2) Add new header name to fstr_ctrl_header_names !
8 ! 3) Describe new function to get parameters from control file !
9 ! in fstr_ctrl.f90 !
10 ! 4) Describe new subroutine to set values of the parameter !
11 ! in fstr_setup.f90 !
12 ! 5) If initial values are necessary, set the value !
13 ! in subroutine fstr_setup_init in fstr_setup.f90 !
15 module m_fstr
16  use hecmw
17  use m_common_struct
18  use m_step
19  use m_out
20  use m_timepoint
21  use mmechgauss
22  use mcontactdef
23 
24  implicit none
25 
26  public
27 
30  integer(kind=kint), parameter :: kyes = 1
31  integer(kind=kint), parameter :: kno = 0
32  integer(kind=kint), parameter :: kon = 1
33  integer(kind=kint), parameter :: koff = 0
34 
36  integer(kind=kint), parameter :: kstprecheck = 0
37  integer(kind=kint), parameter :: kststatic = 1
38  integer(kind=kint), parameter :: ksteigen = 2
39  integer(kind=kint), parameter :: kstheat = 3
40  integer(kind=kint), parameter :: kstdynamic = 4
41  !integer(kind=kint), parameter :: kstNLSTATIC = 5
42  integer(kind=kint), parameter :: kststaticeigen = 6
43  integer(kind=kint), parameter :: kstnzprof = 7
44 
46  integer(kind=kint), parameter :: ksmcg = 1
47  integer(kind=kint), parameter :: ksmbicgstab = 2
48  integer(kind=kint), parameter :: ksmgmres = 3
49  integer(kind=kint), parameter :: ksmgpbicg = 4
50  integer(kind=kint), parameter :: ksmgmresr = 5
51  integer(kind=kint), parameter :: ksmgmresren = 6
52  integer(kind=kint), parameter :: ksmdirect = 101
53 
55  integer(kind=kint), parameter :: knsmnewton = 1
56  integer(kind=kint), parameter :: knsmquasinewton = 2
57 
59  integer(kind=kint), parameter :: kcaslagrange = 1
60  integer(kind=kint), parameter :: kcaalagrange = 2
61 
63  integer(kind=kint), parameter :: kbcffstr = 0 ! BC described in fstr control file (default)
64  integer(kind=kint), parameter :: kbcfnastran = 1 ! nastran file
65 
66  integer(kind=kint), parameter :: kbcinitial = 1
67  integer(kind=kint), parameter :: kbctransit = 2
68 
70  integer(kind=kint), parameter :: restart_outlast = 1
71  integer(kind=kint), parameter :: restart_outall = 2
72 
74  integer(kind=kint), parameter :: kel341fi = 1
75  integer(kind=kint), parameter :: kel341sesns = 2
76 
77  integer(kind=kint), parameter :: kel361fi = 1
78  integer(kind=kint), parameter :: kel361bbar = 2
79  integer(kind=kint), parameter :: kel361ic = 3
80  integer(kind=kint), parameter :: kel361fbar = 4
81 
82  integer(kind=kint), parameter :: kfloadtype_node = 1
83  integer(kind=kint), parameter :: kfloadtype_surf = 2
84 
85  integer(kind=kint), parameter :: kfloadcase_re = 1
86  integer(kind=kint), parameter :: kfloadcase_im = 2
87 
89  integer(kind=kint), parameter :: kitrcontinue = 0
90  integer(kind=kint), parameter :: kitrconverged = 1
91  integer(kind=kint), parameter :: kitrdiverged = 2
92  integer(kind=kint), parameter :: kitrfloatingerror = 3
93 
94 
96  integer(kind = kint) :: myrank
97  integer(kind = kint) :: nprocs
98 
100  logical :: paracontactflag = .false.
101 
103  character(len=HECMW_FILENAME_LEN) :: cntfilname
104  character(len=HECMW_FILENAME_LEN) :: restartfilname
105 
107  integer(kind=kint), parameter :: ilog = 16 ! log
108  integer(kind=kint), parameter :: ista = 17 ! status
109  integer(kind=kint), parameter :: iutb = 18 ! utable
110  integer(kind=kint), parameter :: imsg = 51 ! message (myrank == 0 only)
111  integer(kind=kint), parameter :: idbg = 52 ! debug
112  integer(kind=kint), parameter :: ifvs = 53 ! visual.ini file
113  integer(kind=kint), parameter :: ineu = 54 ! neutral file (heat)
114  integer(kind=kint), parameter :: iresout = 100 ! ~110, keeping for result output file
115 
117  integer(kind=kint) :: sviarray(100)
118  real(kind=kreal) :: svrarray(100)
119 
121  integer(kind=kint), pointer :: iecho
122  integer(kind=kint), pointer :: iresult
123  integer(kind=kint), pointer :: ivisual
124  integer(kind=kint), pointer :: ineutral ! flag for femap neutral file
125  integer(kind=kint), pointer :: irres ! flag for restart, read
126  integer(kind=kint), pointer :: iwres ! flag for restart, write
127  integer(kind=kint), pointer :: nrres ! position of restart read
128  integer(kind=kint), pointer :: nprint ! interval of write
129 
130  integer(kind=kint), parameter :: kopss_solution = 1
131  integer(kind=kint), parameter :: kopss_material = 2
132  integer(kind=kint) :: opsstype = kopss_solution ! output stress/strain type
133 
134 
136  real(kind=kreal), pointer :: ref_temp
137 
139  real(kind=kreal) :: dt ! /=fstr_param%dtime
140  real(kind=kreal) :: etime ! /=fstr_param%etime
141  integer(kind=kint) :: itmax
142  real(kind=kreal) :: eps ! /=fstr_param%eps
143 
145  character(len=HECMW_FILENAME_LEN) :: cond_name
146  integer :: node_elem
147  integer :: grpid
148  integer, pointer :: intval(:) => null()
149  real(kind=kreal), pointer :: realval(:) => null()
150  end type
151  type( tinitialcondition ), pointer, save :: g_initialcnd(:) => null()
152 
155  integer(kind=kint) :: solution_type
156  integer(kind=kint) :: solver_method
157  integer(kind=kint) :: nlsolver_method
158  logical :: nlgeom
159 
161  integer(kind=kint) :: analysis_n
162  real(kind=kreal), pointer :: dtime(:)
163  real(kind=kreal), pointer :: etime(:)
164  real(kind=kreal), pointer :: dtmin(:)
165  real(kind=kreal), pointer :: delmax(:)
166  integer(kind=kint), pointer:: itmax(:)
167  real(kind=kreal), pointer :: eps(:)
168  real(kind=kreal) :: ref_temp
169  integer(kind=kint) :: timepoint_id
170 
172  integer(kind=kint) :: fg_echo
173  integer(kind=kint) :: fg_result
174  integer(kind=kint) :: fg_visual
175 
177  integer(kind=kint) :: fg_neutral
178  integer(kind=kint) :: fg_irres
179  integer(kind=kint) :: fg_iwres
180  integer(kind=kint) :: nrres
181  integer(kind=kint) :: nprint
182 
184  integer(kind=kint) :: n_node
185  integer(kind=kint) :: nn_internal
186  integer(kind=kint), pointer :: global_local_id(:,:)
187 
189  integer( kind=kint ) :: fg_couple
190  integer( kind=kint ) :: fg_couple_type
191  integer( kind=kint ) :: fg_couple_first
192  integer( kind=kint ) :: fg_couple_window
193 
195  integer( kind=kint ) :: restart_out_type
196  integer( kind=kint ) :: restart_version
197 
199  integer( kind=kint ) :: contact_algo
200  type(tcontactparam), pointer :: contactparam(:)
201  type(tcontactinterference), pointer :: contact_if(:)
202 
204  type(tparamautoinc), pointer :: ainc(:)
205  type(time_points), pointer :: timepoints(:)
206  end type fstr_param
207 
209  type( fstr_param ),target :: fstrpr
210 
213  real(kind=kreal), pointer :: stress(:) => null()
214  real(kind=kreal), pointer :: strain(:) => null()
215  real(kind=kreal), pointer :: mises(:) => null()
216 
217  real(kind=kreal), pointer :: pstress(:) => null()
218  real(kind=kreal), pointer :: pstrain(:) => null()
219  real(kind=kreal), pointer :: pstress_vect(:,:) => null()
220  real(kind=kreal), pointer :: pstrain_vect(:,:) => null()
221 
222  real(kind=kreal), pointer :: estress(:) => null()
223  real(kind=kreal), pointer :: estrain(:) => null()
224  real(kind=kreal), pointer :: emises(:) => null()
225  real(kind=kreal), pointer :: eplstrain(:) => null()
226 
227  real(kind=kreal), pointer :: epstress(:) => null()
228  real(kind=kreal), pointer :: epstrain(:) => null()
229  real(kind=kreal), pointer :: epstress_vect(:,:) => null()
230  real(kind=kreal), pointer :: epstrain_vect(:,:) => null()
231  real(kind=kreal), pointer :: enqm(:) => null()
232 
233 
234  type(fstr_solid_physic_val), pointer :: layer(:) => null()
235  type(fstr_solid_physic_val), pointer :: plus => null()
236  type(fstr_solid_physic_val), pointer :: minus => null()
237  end type fstr_solid_physic_val
238 
240  integer(kind=kint) :: file_type ! kbcfFSTR or kbcfNASTRAN
241  integer(kind=kint) :: statictype ! 1:Total, 2:Updated, 3:Infinitesimal
242  integer(kind=kint) :: nstep_tot
243 
244  type(step_info), pointer :: step_ctrl(:) =>null()
245  type(t_output_ctrl), pointer :: output_ctrl(:)=>null()
246 
248  integer(kind=kint) :: boundary_ngrp_tot
249  integer(kind=kint), pointer :: boundary_ngrp_grpid (:) =>null()
250  integer(kind=kint), pointer :: boundary_ngrp_id (:) =>null()
251  integer(kind=kint), pointer :: boundary_ngrp_type (:) =>null()
252  integer(kind=kint), pointer :: boundary_ngrp_amp (:) =>null()
253  real(kind=kreal), pointer :: boundary_ngrp_val(:) =>null()
254  integer(kind=kint), pointer :: boundary_ngrp_istot (:) =>null()
255  integer(kind=kint) :: boundary_ngrp_rot
256  integer(kind=kint), pointer :: boundary_ngrp_rotid (:) =>null()
257  integer(kind=kint), pointer :: boundary_ngrp_centerid (:) =>null()
258 
260  integer(kind=kint) :: velocity_type
261  integer(kind=kint) :: velocity_ngrp_tot
262  integer(kind=kint), pointer :: velocity_ngrp_grpid (:) =>null()
263  integer(kind=kint), pointer :: velocity_ngrp_id (:) =>null()
264  integer(kind=kint), pointer :: velocity_ngrp_type (:) =>null()
265  integer(kind=kint), pointer :: velocity_ngrp_amp (:) =>null()
266  real(kind=kreal), pointer :: velocity_ngrp_val(:) =>null()
267 
269  integer(kind=kint) :: acceleration_type
270  integer(kind=kint) :: acceleration_ngrp_tot
271  integer(kind=kint), pointer :: acceleration_ngrp_grpid (:) =>null()
272  integer(kind=kint), pointer :: acceleration_ngrp_id (:) =>null()
273  integer(kind=kint), pointer :: acceleration_ngrp_type (:) =>null()
274  integer(kind=kint), pointer :: acceleration_ngrp_amp (:) =>null()
275  real(kind=kreal), pointer :: acceleration_ngrp_val(:) =>null()
276 
278  integer(kind=kint) :: cload_ngrp_tot
279  integer(kind=kint), pointer :: cload_ngrp_grpid (:) =>null()
280  integer(kind=kint), pointer :: cload_ngrp_id (:)
281  integer(kind=kint), pointer :: cload_ngrp_dof (:)
282  integer(kind=kint), pointer :: cload_ngrp_amp (:)
283  real(kind=kreal), pointer :: cload_ngrp_val(:)
284  integer(kind=kint) :: cload_ngrp_rot
285  integer(kind=kint), pointer :: cload_ngrp_rotid (:) =>null()
286  integer(kind=kint), pointer :: cload_ngrp_centerid (:) =>null()
287 
289  integer(kind=kint) :: dload_ngrp_tot
290  integer(kind=kint) :: dload_follow
291  integer(kind=kint), pointer :: dload_ngrp_grpid (:) =>null()
292  integer(kind=kint), pointer :: dload_ngrp_id (:)
293  integer(kind=kint), pointer :: dload_ngrp_lid (:)
294  integer(kind=kint), pointer :: dload_ngrp_amp (:)
295  real(kind=kreal), pointer :: dload_ngrp_params(:,:)
296 
298  integer(kind=kint) :: temp_ngrp_tot
299  integer(kind=kint) :: temp_irres
300  integer(kind=kint) :: temp_tstep
301  integer(kind=kint) :: temp_interval
302  integer(kind=kint) :: temp_rtype ! type of reading result; 1: step-based; 2: time-based
303  real(kind=kreal) :: temp_factor
304  integer(kind=kint), pointer :: temp_ngrp_grpid (:) =>null()
305  integer(kind=kint), pointer :: temp_ngrp_id (:)
306  real(kind=kreal), pointer :: temp_ngrp_val(:)
307 
309  integer(kind=kint) :: spring_ngrp_tot
310  integer(kind=kint), pointer :: spring_ngrp_grpid (:) =>null()
311  integer(kind=kint), pointer :: spring_ngrp_id (:)
312  integer(kind=kint), pointer :: spring_ngrp_dof (:)
313  integer(kind=kint), pointer :: spring_ngrp_amp (:)
314  real(kind=kreal), pointer :: spring_ngrp_val(:)
315 
317  integer( kind=kint ) :: couple_ngrp_tot
318  integer( kind=kint ),pointer :: couple_ngrp_id(:)
319 
321  integer(kind=kint) :: maxn_gauss
322 
323  real(kind=kreal), pointer :: stress(:)
324  real(kind=kreal), pointer :: strain(:)
325  real(kind=kreal), pointer :: mises(:)
326 
327  real(kind=kreal), pointer :: pstress(:)
328  real(kind=kreal), pointer :: pstrain(:)
329  real(kind=kreal), pointer :: pstress_vect(:,:)
330  real(kind=kreal), pointer :: pstrain_vect(:,:)
331 
332  real(kind=kreal), pointer :: estress(:)
333  real(kind=kreal), pointer :: estrain(:)
334  real(kind=kreal), pointer :: emises(:)
335 
336  real(kind=kreal), pointer :: epstress(:)
337  real(kind=kreal), pointer :: epstrain(:)
338  real(kind=kreal), pointer :: epstress_vect(:,:)
339  real(kind=kreal), pointer :: epstrain_vect(:,:)
340 
341  real(kind=kreal), pointer :: tnstrain(:)
342  real(kind=kreal), pointer :: testrain(:)
343 
344  real(kind=kreal), pointer :: yield_ratio(:)
345 
346  real(kind=kreal), pointer :: enqm(:)
347  real(kind=kreal), pointer :: reaction(:)
348 
349  real(kind=kreal), pointer :: cont_nforce(:)
350  real(kind=kreal), pointer :: cont_fric(:)
351  real(kind=kreal), pointer :: cont_relvel(:)
352  real(kind=kreal), pointer :: cont_state(:)
353  integer(kind=kint), pointer :: cont_sgrp_id(:)
354  real(kind=kreal), pointer :: cont_area(:)
355  real(kind=kreal), pointer :: cont_ntrac(:)
356  real(kind=kreal), pointer :: cont_ftrac(:)
357  real(kind=kreal), pointer :: embed_nforce(:)
358 
359  type(fstr_solid_physic_val), pointer :: solid=>null()
360  type(fstr_solid_physic_val), pointer :: shell=>null()
361  type(fstr_solid_physic_val), pointer :: beam =>null()
362 
364  integer(kind=kint) :: restart_nout
367  integer(kind=kint) :: restart_nin
368  type(step_info) :: step_ctrl_restart
369 
370  integer(kind=kint) :: max_lyr
371  integer(kind=kint) :: is_33shell
372  integer(kind=kint) :: is_33beam
373  integer(kind=kint) :: is_heat
374  integer(kind=kint) :: max_ncon_stf
375  integer(kind=kint) :: max_ncon
376  integer(kind=kint), pointer :: is_rot(:) => null()
377  integer(kind=kint) :: elemopt361
378  logical :: is_smoothing_active
379  real(kind=kreal) :: factor(2)
382  integer(kind=kint) :: nrstat_i(10)
383  real(kind=kreal) :: nrstat_r(10)
384  integer(kind=kint) :: autoinc_stat
385  integer(kind=kint) :: cutback_stat
386 
387  real(kind=kreal), pointer :: gl(:)
388  real(kind=kreal), pointer :: gl0(:)
389  real(kind=kreal), pointer :: eforce(:)
390  real(kind=kreal), pointer :: qforce(:)
391  real(kind=kreal), pointer :: qforce_bak(:)
392  real(kind=kreal), pointer :: unode(:) => null()
393  real(kind=kreal), pointer :: unode_bak(:) => null()
394  real(kind=kreal), pointer :: dunode(:) => null()
395  real(kind=kreal), pointer :: ddunode(:) => null()
396  real(kind=kreal), pointer :: temperature(:)=> null()
397  real(kind=kreal), pointer :: temp_bak(:) => null()
398  real(kind=kreal), pointer :: last_temp(:) => null()
399 
400  type( telement ), pointer :: elements(:) =>null()
401  type( tmaterial ),pointer :: materials(:) =>null()
402  integer :: n_contacts
403  type( tcontact ), pointer :: contacts(:) =>null()
404  integer :: n_embeds
405  type( tcontact ), pointer :: embeds(:) =>null()
406  integer :: n_fix_mpc
407  real(kind=kreal), pointer :: mpc_const(:) =>null()
408  type(tsection), pointer :: sections(:) =>null()
409 
410  ! for cutback
411  ! ####################### Notice #######################
412  ! # If you add new variables to store analysis status, #
413  ! # - backup variables with postfix "_bkup" here #
414  ! # - backup process to module m_fstr_Cutback #
415  ! # must be added if necessary. #
416  ! ######################################################
417  real(kind=kreal), pointer :: unode_bkup(:) => null()
418  real(kind=kreal), pointer :: qforce_bkup(:) => null()
419  real(kind=kreal), pointer :: last_temp_bkup(:) => null()
420  type( telement ), pointer :: elements_bkup(:) =>null()
421  type( tcontact ), pointer :: contacts_bkup(:) =>null()
422  type( tcontact ), pointer :: embeds_bkup(:) =>null()
423  end type fstr_solid
424 
428  integer(kind=kint) :: is_steady
429  real(kind=kreal) :: beta
430  logical :: is_iter_max_limit
431 
433  integer(kind=kint) :: steptot
434  integer(kind=kint) :: restart_nout
435  real(kind=kreal), pointer :: step_dltime(:), step_eetime(:)
436  real(kind=kreal), pointer :: step_delmin(:), step_delmax(:)
437  integer(kind=kint) :: timepoint_id
438 
440  integer(kind=kint) :: materialtot
441  integer(kind=kint), pointer :: rhotab(:), cptab(:), condtab(:)
442  real(kind=kreal), pointer :: rho(:,:), rhotemp(:,:)
443  real(kind=kreal), pointer :: cp(:,:), cptemp(:,:)
444  real(kind=kreal), pointer :: cond(:,:),condtemp(:,:)
445 
446  real(kind=kreal), pointer :: rhofunca(:,:), rhofuncb(:,:)
447  real(kind=kreal), pointer :: cpfunca(:,:), cpfuncb(:,:)
448  real(kind=kreal), pointer :: condfunca(:,:),condfuncb(:,:)
449 
451  integer(kind=kint) :: amplitudetot
452  integer(kind=kint), pointer :: ampltab(:)
453  real(kind=kreal), pointer :: ampl(:,:), ampltime(:,:)
454  real(kind=kreal), pointer :: amplfunca(:,:), amplfuncb(:,:)
455 
457  real(kind=kreal), pointer :: temp0(:)
458  real(kind=kreal), pointer :: tempc(:)
459  real(kind=kreal), pointer :: temp(:)
460 
462  integer(kind=kint) :: t_fix_tot
463  integer(kind=kint), pointer :: t_fix_node(:)
464  integer(kind=kint), pointer :: t_fix_ampl(:)
465  real(kind=kreal), pointer :: t_fix_val(:)
466 
468  integer(kind=kint) :: q_nod_tot
469  integer(kind=kint), pointer :: q_nod_node(:)
470  integer(kind=kint), pointer :: q_nod_ampl(:)
471  real(kind=kreal), pointer :: q_nod_val(:)
472 
474  integer(kind=kint) :: q_vol_tot
475  integer(kind=kint), pointer :: q_vol_elem(:)
476  integer(kind=kint), pointer :: q_vol_ampl(:)
477  real(kind=kreal), pointer :: q_vol_val(:)
478 
480  integer(kind=kint) :: q_suf_tot
481  integer(kind=kint), pointer :: q_suf_elem(:)
482  integer(kind=kint), pointer :: q_suf_ampl(:)
483  integer(kind=kint), pointer :: q_suf_surf(:)
484  real(kind=kreal), pointer :: q_suf_val(:)
485 
487  integer(kind=kint) :: r_suf_tot
488  integer(kind=kint), pointer :: r_suf_elem(:)
489  integer(kind=kint), pointer :: r_suf_ampl(:,:)
490  integer(kind=kint), pointer :: r_suf_surf(:)
491  real(kind=kreal), pointer :: r_suf_val(:,:)
492 
494  integer(kind=kint) :: h_suf_tot
495  integer(kind=kint), pointer :: h_suf_elem(:)
496  integer(kind=kint), pointer :: h_suf_ampl(:,:)
497  integer(kind=kint), pointer :: h_suf_surf(:)
498  real(kind=kreal), pointer :: h_suf_val(:,:)
499 
500  integer(kind=kint) :: wl_tot
501  type(tweldline), pointer :: weldline(:) => null()
502  end type fstr_heat
503 
507  integer(kind=kint) :: idx_eqa ! implicit or explicit
508  integer(kind=kint) :: idx_resp ! time history or steady-state harmonic response analysis
509 
511  integer(kind=kint) :: n_step ! total step number of analysis
512  real(kind=kreal) :: t_start ! start time of analysis
513  real(kind=kreal) :: t_curr ! current time of analysis
514  real(kind=kreal) :: t_end ! end time of analysis
515  real(kind=kreal) :: t_delta ! time increment
516  integer(kind=kint) :: restart_nout ! output interval of restart file
517  ! (if .gt.0) restart file write
518  ! (if .lt.0) restart file read and write
519  integer(kind=kint) :: restart_nin !input number of restart file
520 
522  real(kind=kreal) :: ganma ! Newmark-beta parameter ganma
523  real(kind=kreal) :: beta ! Newmark-beta parameter beta
524 
526  integer(kind=kint) :: idx_mas ! mass matrix type
527 
529  integer(kind=kint) :: idx_dmp ! damping type
530  real(kind=kreal) :: ray_m ! Rayleigh damping parameter Rm
531  real(kind=kreal) :: ray_k ! Rayleigh damping parameter Rk
532 
534  logical :: varinitialize ! initialization flag
535 
537  integer(kind=kint) :: nout ! output interval of result
538  integer(kind=kint) :: ngrp_monit ! node of monitoring result
539  integer(kind=kint) :: nout_monit ! output interval of result monitoring
540  integer(kind=kint) :: i_step ! step number
541  integer(kind=kint) :: iout_list(6) ! 0:not output 1:output
542  ! iout_list(1): displacement
543  ! iout_list(2): velocity
544  ! iout_list(3): acceleration
545  ! iout_list(4): reaction force
546  ! iout_list(5): strain
547  ! iout_list(6): stress
548 
550  real(kind=kreal), pointer :: disp(:,:)
551  real(kind=kreal), pointer :: vel(:,:)
552  real(kind=kreal), pointer :: acc(:,:)
553 
554  real(kind=kreal) :: kineticenergy
555  real(kind=kreal) :: strainenergy
556  real(kind=kreal) :: totalenergy
557 
559  real(kind=kreal), pointer :: vec1(:)
560  real(kind=kreal), pointer :: vec2(:)
561  real(kind=kreal), pointer :: vec3(:)
562 
563  integer(kind=kint) :: dynamic_iw4 = 204
564  integer(kind=kint) :: dynamic_iw5 = 205
565  integer(kind=kint) :: dynamic_iw6 = 206
566  integer(kind=kint) :: dynamic_iw7 = 207
567  integer(kind=kint) :: dynamic_iw8 = 208
568  integer(kind=kint) :: dynamic_iw9 = 209
569  integer(kind=kint) :: dynamic_iw10 = 210
570  end type fstr_dynamic
571 
573  integer(kind=kint) :: fload_ngrp_tot
574  integer(kind=kint), pointer :: fload_ngrp_grpid(:) => null()
575  integer(kind=kint), pointer :: fload_ngrp_id(:) => null()
576  integer(kind=kint), pointer :: fload_ngrp_type(:) => null()
577  integer(kind=kint), pointer :: fload_ngrp_dof(:) => null()
578  real(kind=kreal), pointer :: fload_ngrp_valre(:) => null()
579  real(kind=kreal), pointer :: fload_ngrp_valim(:) => null()
580  character(len=HECMW_FILENAME_LEN) :: eigenlog_filename
581  integer(kind=kint) :: start_mode
582  integer(kind=kint) :: end_mode
583  end type fstr_freqanalysis
584 
586  integer(kind=kint) :: nummode
587  integer(kind=kint) :: numnodedof
588  real(kind=kreal), pointer :: eigomega(:) => null()
589  real(kind=kreal), pointer :: eigvector(:,:) => null()
590  real(kind=kreal) :: rayalpha, raybeta
591  end type fstr_freqanalysis_data
592 
596  integer(kind=kint) :: nget ! Solved eigen value number (default:5)
597  integer(kind=kint) :: maxiter ! Max. Lcz iterations (default:60)
598  integer(kind=kint) :: iter ! Max. Lcz iterations (default:60)
599  real (kind=kreal) :: sigma ! 0.0
600  real (kind=kreal) :: tolerance ! Lcz tolerance (default:1.0e-8)
601  real (kind=kreal) :: totalmass
602  real (kind=kreal), pointer :: eigval(:)
603  real (kind=kreal), pointer :: eigvec(:,:)
604  real (kind=kreal), pointer :: filter(:)
605  real (kind=kreal), pointer :: mass(:)
606  real (kind=kreal), pointer :: effmass(:)
607  real (kind=kreal), pointer :: partfactor(:)
608  logical :: is_free = .false.
609  end type fstr_eigen
610 
614  integer( kind=kint ) :: dof ! == 3
615  integer( kind=kint ) :: ndof ! total dof (coupled_node_n*dof)
616  integer( kind=kint ) :: coupled_node_n
618  integer, pointer :: coupled_node(:) ! local node id sent from revocap
619  real( kind=8 ),pointer :: trac(:) ! input (x,y,z,x,y,z ... )
620  real( kind=8 ),pointer :: disp(:) ! output (x,y,z,x,y,z ... )
621  real( kind=8 ),pointer :: velo(:) ! output (x,y,z,x,y,z ... )
622  real( kind=8 ),pointer :: accel(:) ! output (x,y,z,x,y,z ... )
624  integer( kind=kint ),pointer :: index(:) ! size:total node num.
626  end type fstr_couple
627 
630  integer(kind=kint) :: egrpid
631  real( kind=kreal ) :: i
632  real( kind=kreal ) :: u
633  real( kind=kreal ) :: coe
634  real( kind=kreal ) :: v
635  integer(kind=kint) :: xyz
636  real(kind=kreal) :: n1, n2
637  real(kind=kreal) :: distol
638  real(kind=kreal) :: tstart
639  end type tweldline
640 
642  type tsection
643  !integer :: mat_ID
644  !integer :: iset
645  !integer :: orien_ID
646  !real(kind=kreal) :: thickness
647  integer :: elemopt341
648  !integer :: elemopt342
649  !integer :: elemopt351
650  !integer :: elemopt352
651  integer :: elemopt361
652  !integer :: elemopt362
653  end type tsection
654 
655 contains
656 
658  subroutine fstr_nullify_fstr_param( P )
659  implicit none
660  type( fstr_param ) :: P
661 
662  nullify( p%dtime )
663  nullify( p%etime )
664  nullify( p%dtmin )
665  nullify( p%delmax )
666  nullify( p%itmax )
667  nullify( p%eps )
668  nullify( p%global_local_ID)
669  nullify( p%timepoints )
670  end subroutine fstr_nullify_fstr_param
671 
672  subroutine fstr_nullify_fstr_solid( S )
673  implicit none
674  type( fstr_solid ) :: S
675 
676  nullify( s%BOUNDARY_ngrp_ID )
677  nullify( s%BOUNDARY_ngrp_type )
678  nullify( s%BOUNDARY_ngrp_amp )
679  nullify( s%BOUNDARY_ngrp_val)
680  nullify( s%BOUNDARY_ngrp_rotID )
681  nullify( s%BOUNDARY_ngrp_centerID )
682  nullify( s%CLOAD_ngrp_ID )
683  nullify( s%CLOAD_ngrp_DOF )
684  nullify( s%CLOAD_ngrp_amp )
685  nullify( s%CLOAD_ngrp_rotID )
686  nullify( s%CLOAD_ngrp_centerID )
687  nullify( s%CLOAD_ngrp_val )
688  nullify( s%DLOAD_ngrp_ID )
689  nullify( s%DLOAD_ngrp_LID )
690  nullify( s%DLOAD_ngrp_amp )
691  nullify( s%DLOAD_ngrp_params )
692  nullify( s%TEMP_ngrp_ID )
693  nullify( s%TEMP_ngrp_val )
694  nullify( s%SPRING_ngrp_ID )
695  nullify( s%SPRING_ngrp_DOF )
696  nullify( s%SPRING_ngrp_amp )
697  nullify( s%SPRING_ngrp_val )
698  nullify( s%STRESS )
699  nullify( s%STRAIN )
700  nullify( s%MISES )
701  nullify( s%PSTRESS )
702  nullify( s%PSTRAIN )
703  nullify( s%PSTRESS_VECT )
704  nullify( s%PSTRAIN_VECT )
705  nullify( s%REACTION )
706  nullify( s%ESTRESS )
707  nullify( s%ESTRAIN )
708  nullify( s%EMISES )
709  nullify( s%EPSTRESS )
710  nullify( s%EPSTRAIN )
711  nullify( s%EPSTRESS_VECT )
712  nullify( s%EPSTRAIN_VECT )
713  nullify( s%ENQM )
714  nullify( s%GL )
715  nullify( s%GL0 )
716  nullify( s%QFORCE )
717  nullify( s%VELOCITY_ngrp_ID )
718  nullify( s%VELOCITY_ngrp_type )
719  nullify( s%VELOCITY_ngrp_amp )
720  nullify( s%VELOCITY_ngrp_val )
721  nullify( s%ACCELERATION_ngrp_ID )
722  nullify( s%ACCELERATION_ngrp_type )
723  nullify( s%ACCELERATION_ngrp_amp )
724  nullify( s%ACCELERATION_ngrp_val )
725  nullify( s%COUPLE_ngrp_ID )
726  end subroutine fstr_nullify_fstr_solid
727 
728  subroutine fstr_nullify_fstr_heat( H )
729  implicit none
730  type( fstr_heat ) :: H
731 
732  nullify( h%STEP_DLTIME )
733  nullify( h%STEP_EETIME )
734  nullify( h%STEP_DELMIN )
735  nullify( h%STEP_DELMAX )
736  nullify( h%RHO )
737  nullify( h%RHOtemp )
738  nullify( h%CP )
739  nullify( h%CPtemp )
740  nullify( h%COND )
741  nullify( h%CONDtemp )
742  nullify( h%RHOtab )
743  nullify( h%CPtab )
744  nullify( h%CONDtab )
745  nullify( h%RHOfuncA )
746  nullify( h%RHOfuncB )
747  nullify( h%CPfuncA )
748  nullify( h%CPfuncB )
749  nullify( h%CONDfuncA )
750  nullify( h%CONDfuncB )
751  nullify( h%AMPL )
752  nullify( h%AMPLtime )
753  nullify( h%AMPLtab )
754  nullify( h%AMPLfuncA )
755  nullify( h%AMPLfuncB )
756  nullify( h%TEMP0 )
757  nullify( h%TEMPC )
758  nullify( h%TEMP )
759  nullify( h%T_FIX_node )
760  nullify( h%T_FIX_ampl )
761  nullify( h%T_FIX_val )
762  nullify( h%Q_NOD_node )
763  nullify( h%Q_NOD_ampl )
764  nullify( h%Q_NOD_val )
765  nullify( h%Q_VOL_elem )
766  nullify( h%Q_VOL_ampl )
767  nullify( h%Q_VOL_val )
768  nullify( h%Q_SUF_elem )
769  nullify( h%Q_SUF_ampl )
770  nullify( h%Q_SUF_surf )
771  nullify( h%Q_SUF_val )
772  nullify( h%R_SUF_elem )
773  nullify( h%R_SUF_ampl )
774  nullify( h%R_SUF_surf )
775  nullify( h%R_SUF_val )
776  nullify( h%H_SUF_elem )
777  nullify( h%H_SUF_ampl )
778  nullify( h%H_SUF_surf )
779  nullify( h%H_SUF_val )
780  end subroutine fstr_nullify_fstr_heat
781 
782  subroutine fstr_nullify_fstr_dynamic( DY )
783  implicit none
784  type( fstr_dynamic ) :: DY
785 
786  nullify( dy%DISP )
787  nullify( dy%VEL )
788  nullify( dy%ACC )
789  nullify( dy%VEC1 )
790  nullify( dy%VEC2 )
791  nullify( dy%VEC3 )
792  end subroutine fstr_nullify_fstr_dynamic
793 
794  subroutine fstr_nullify_fstr_freqanalysis( f )
795  implicit none
796  type( fstr_freqanalysis ), intent(inout) :: f
797 
798  f%FLOAD_ngrp_tot = 0
799  nullify( f%FLOAD_ngrp_GRPID )
800  nullify( f%FLOAD_ngrp_ID )
801  nullify( f%FLOAD_ngrp_TYPE )
802  nullify( f%FLOAD_ngrp_DOF )
803  nullify( f%FLOAD_ngrp_valre )
804  nullify( f%FLOAD_ngrp_valim )
805  end subroutine fstr_nullify_fstr_freqanalysis
806 
807  subroutine fstr_nullify_fstr_eigen( E )
808  implicit none
809  type( fstr_eigen ) :: E
810 
811  nullify( e%mass )
812  end subroutine fstr_nullify_fstr_eigen
813 
814  subroutine fstr_nullify_fstr_couple( C )
815  implicit none
816  type( fstr_couple ) :: C
817 
818  nullify( c%coupled_node )
819  nullify( c%trac )
820  nullify( c%velo )
821  nullify( c%index )
822  end subroutine fstr_nullify_fstr_couple
823 
825  subroutine fstr_mat_init( hecMAT )
826  implicit none
827  type(hecmwst_matrix) :: hecMAT
828 
829  hecmat%Iarray(1) = 100 ! = nier
830  hecmat%Iarray(2) = 1 ! = method
831  hecmat%Iarray(3) = 1 ! = precond
832  hecmat%Iarray(4) = 0 ! = nset
833  hecmat%Iarray(5) = 1 ! = iterpremax
834  hecmat%Iarray(6) = 10 ! = nrest
835  hecmat%Iarray(7) = 0 ! = scaling
836  hecmat%Iarray(21)= kno ! = iterlog
837  hecmat%Iarray(22)= kno ! = timelog
838  hecmat%Iarray(31)= 0 ! = dumptype
839  hecmat%Iarray(32)= 0 ! = dumpexit
840  hecmat%Iarray(33)= 0 ! = usejad
841  hecmat%Iarray(34)= 10 ! = ncolor_in
842  hecmat%Iarray(13)= 0 ! = mpc_method
843  hecmat%Iarray(14)= 0 ! = estcond
844  hecmat%Iarray(35)= 3 ! = maxrecycle_precond
845  hecmat%Iarray(41)= 0 ! = solver_opt1
846  hecmat%Iarray(42)= 0 ! = solver_opt2
847  hecmat%Iarray(43)= 0 ! = solver_opt3
848  hecmat%Iarray(44)= 0 ! = solver_opt4
849  hecmat%Iarray(45)= 0 ! = solver_opt5
850  hecmat%Iarray(46)= 0 ! = solver_opt6
851 
852  hecmat%Rarray(1) = 1.0e-8 ! = resid
853  hecmat%Rarray(2) = 1.0 ! = sigma_diag
854  hecmat%Rarray(3) = 0.0 ! = sigma
855  hecmat%Rarray(4) = 0.1 ! = thresh
856  hecmat%Rarray(5) = 0.1 ! = filter
857  hecmat%Rarray(11)= 1.0e+4 ! = penalty
858 
859  hecmat%Iarray(96) = 0 ! nrecycle_precond
860  hecmat%Iarray(97) = kyes ! flag_numfact
861  hecmat%Iarray(98) = kyes ! flag_symbfact
862  hecmat%Iarray(99) = kyes ! indirect method
863  end subroutine fstr_mat_init
864 
865  subroutine hecmat_init( hecMAT )
866  implicit none
867  type( hecmwst_matrix ) :: hecMAT
868  integer :: ndof, nn, ierror
869  ndof = hecmat%NDOF
870  nn = ndof*ndof
871  allocate (hecmat%AL(nn*hecmat%NPL) ,stat=ierror )
872  if( ierror /= 0 ) then
873  write(*,*) "##ERROR : not enough memory"
874  write(idbg,*) 'stop due to allocation error'
875  call flush(idbg)
876  call hecmw_abort( hecmw_comm_get_comm() )
877  end if
878  allocate (hecmat%AU(nn*hecmat%NPU) ,stat=ierror )
879  if( ierror /= 0 ) then
880  write(*,*) "##ERROR : not enough memory"
881  write(idbg,*) 'stop due to allocation error'
882  call flush(idbg)
883  call hecmw_abort( hecmw_comm_get_comm() )
884  end if
885  allocate (hecmat%B(ndof*hecmat%NP) ,stat=ierror )
886  if( ierror /= 0 ) then
887  write(*,*) "##ERROR : not enough memory"
888  write(idbg,*) 'stop due to allocation error'
889  call flush(idbg)
890  call hecmw_abort( hecmw_comm_get_comm() )
891  end if
892  hecmat%B(:)=0.d0
893  allocate (hecmat%D(nn*hecmat%NP) ,stat=ierror )
894  if( ierror /= 0 ) then
895  write(*,*) "##ERROR : not enough memory"
896  write(idbg,*) 'stop due to allocation error'
897  call flush(idbg)
898  call hecmw_abort( hecmw_comm_get_comm() )
899  end if
900  allocate (hecmat%X(ndof*hecmat%NP) ,stat=ierror )
901  if( ierror /= 0 ) then
902  write(*,*) "##ERROR : not enough memory"
903  write(idbg,*) 'stop due to allocation error'
904  call flush(idbg)
905  call hecmw_abort( hecmw_comm_get_comm() )
906  end if
907  allocate (hecmat%ALU(nn*hecmat%N) ,stat=ierror )
908  if( ierror /= 0 ) then
909  write(*,*) "##ERROR : not enough memory"
910  write(idbg,*) 'stop due to allocation error'
911  call flush(idbg)
912  call hecmw_abort( hecmw_comm_get_comm() )
913  endif
914  hecmat%D = 0.0d0
915  hecmat%AL = 0.0d0
916  hecmat%AU = 0.0d0
917  hecmat%B = 0.0d0
918  hecmat%X = 0.0d0
919  hecmat%ALU = 0.0d0
920  end subroutine hecmat_init
921 
922  subroutine hecmat_finalize( hecMAT )
923  implicit none
924  type( hecmwst_matrix ) :: hecMAT
925  integer :: ndof, nn, ierror
926  ndof = hecmat%NDOF
927  nn = ndof*ndof
928  if( associated(hecmat%AL) ) then
929  deallocate(hecmat%AL ,stat=ierror)
930  if( ierror /= 0 ) then
931  write(idbg,*) 'stop due to deallocation error'
932  call flush(idbg)
933  call hecmw_abort( hecmw_comm_get_comm())
934  end if
935  endif
936  if( associated(hecmat%AU) ) then
937  deallocate(hecmat%AU ,stat=ierror)
938  if( ierror /= 0 ) then
939  write(idbg,*) 'stop due to deallocation error'
940  call flush(idbg)
941  call hecmw_abort( hecmw_comm_get_comm())
942  end if
943  endif
944  if( associated(hecmat%B) ) then
945  deallocate(hecmat%B ,stat=ierror)
946  if( ierror /= 0 ) then
947  write(idbg,*) 'stop due to deallocation error'
948  call flush(idbg)
949  call hecmw_abort( hecmw_comm_get_comm())
950  end if
951  endif
952  if( associated(hecmat%D) ) then
953  deallocate(hecmat%D ,stat=ierror)
954  if( ierror /= 0 ) then
955  write(idbg,*) 'stop due to deallocation error'
956  call flush(idbg)
957  call hecmw_abort( hecmw_comm_get_comm())
958  end if
959  endif
960  if( associated(hecmat%X) ) then
961  deallocate(hecmat%X ,stat=ierror)
962  if( ierror /= 0 ) then
963  write(idbg,*) 'stop due to deallocation error'
964  call flush(idbg)
965  call hecmw_abort( hecmw_comm_get_comm())
966  end if
967  endif
968  if( associated(hecmat%ALU) ) then
969  deallocate(hecmat%ALU ,stat=ierror)
970  if( ierror /= 0 ) then
971  write(idbg,*) 'stop due to deallocation error'
972  call flush(idbg)
973  call hecmw_abort( hecmw_comm_get_comm())
974  end if
975  endif
976  end subroutine hecmat_finalize
977 
979  subroutine fstr_param_init( fstrPARAM, hecMESH )
980  implicit none
981  type(fstr_param) :: fstrPARAM
982  type(hecmwst_local_mesh) :: hecMESH
983  integer(kind=kint) :: i
984  external fstr_sort_index
985 
986  fstrparam%solution_type = kststatic
987  fstrparam%nlgeom = .false.
988  fstrparam%solver_method = ksmcg
989  fstrparam%nlsolver_method = knsmnewton
990 
991  !!STATIC !HEAT
992  fstrparam%analysis_n = 0
993  fstrparam%ref_temp = 0
994 
995  ! output control
996  fstrparam%fg_echo = koff
997  fstrparam%fg_result = koff
998  fstrparam%fg_visual = koff
999 
1000  ! for heat ...
1001  fstrparam%fg_neutral = koff
1002  fstrparam%fg_irres = kno
1003  fstrparam%fg_iwres = kno
1004  fstrparam%nrres = 1
1005  fstrparam%nprint = 1
1006 
1007  ! for couple
1008  fstrparam%fg_couple = 0
1009  fstrparam%fg_couple_type = 0
1010  fstrparam%fg_couple_first= 0
1011  fstrparam%fg_couple_window= 0
1012 
1013  ! for restart control
1014  fstrparam%restart_version = 5
1015 
1016  ! index table for global node ID sorting
1017 
1018  fstrparam%n_node = hecmesh%n_node;
1019  fstrparam%nn_internal = hecmesh%nn_internal;
1020  allocate( fstrparam%global_local_ID(2,hecmesh%nn_internal))
1021  do i = 1, hecmesh%nn_internal
1022  fstrparam%global_local_ID(1,i) = hecmesh%global_node_ID(i)
1023  fstrparam%global_local_ID(2,i) = i
1024  end do
1025  call fstr_sort_index( fstrparam%global_local_ID, hecmesh%nn_internal )
1026  end subroutine fstr_param_init
1027 
1028  logical function fstr_isboundaryactive( fstrSOLID, nbc, cstep )
1029  type(fstr_solid) :: fstrsolid
1030  integer, intent(in) :: nbc
1031  integer, intent(in) :: cstep
1032  fstr_isboundaryactive = .true.
1033  if( .not. associated(fstrsolid%step_ctrl) ) return
1034  if( cstep>fstrsolid%nstep_tot ) return
1035  fstr_isboundaryactive = isboundaryactive( nbc, fstrsolid%step_ctrl(cstep) )
1036  end function
1037 
1038  logical function fstr_isloadactive( fstrSOLID, nbc, cstep )
1039  type(fstr_solid) :: fstrsolid
1040  integer, intent(in) :: nbc
1041  integer, intent(in) :: cstep
1042  fstr_isloadactive = .true.
1043  if( cstep > 0 ) then
1044  if( .not. associated(fstrsolid%step_ctrl) ) return
1045  if( cstep>fstrsolid%nstep_tot ) return
1046  fstr_isloadactive = isloadactive( nbc, fstrsolid%step_ctrl(cstep) )
1047  else
1048  fstr_isloadactive = isloadactive( nbc, fstrsolid%step_ctrl_restart )
1049  endif
1050  end function
1051 
1052  logical function fstr_iscontactactive( fstrSOLID, nbc, cstep )
1053  type(fstr_solid) :: fstrsolid
1054  integer, intent(in) :: nbc
1055  integer, intent(in) :: cstep
1056  fstr_iscontactactive = .true.
1057  if( .not. associated(fstrsolid%step_ctrl) ) return
1058  if( cstep>fstrsolid%nstep_tot ) return
1059  fstr_iscontactactive = iscontactactive( nbc, fstrsolid%step_ctrl(cstep) )
1060  end function
1061 
1062  logical function fstr_isembedactive( fstrSOLID, nbc, cstep )
1063  type(fstr_solid) :: fstrsolid
1064  integer, intent(in) :: nbc
1065  integer, intent(in) :: cstep
1066  fstr_isembedactive = .true.
1067  if( .not. associated(fstrsolid%step_ctrl) ) return
1068  if( cstep>fstrsolid%nstep_tot ) return
1069  fstr_isembedactive = iscontactactive( nbc, fstrsolid%step_ctrl(cstep) )
1070  end function
1071 
1073  subroutine get_coordsys( cdsys_ID, hecMESH, fstrSOLID, coords )
1074  integer, intent(in) :: cdsys_ID
1075  type(hecmwst_local_mesh) :: hecMESH
1076  type(fstr_solid), intent(inout) :: fstrSOLID
1077  real(kind=kreal), intent(out) :: coords(3,3)
1078  integer :: ik
1079 
1080  coords = 0.d0
1081  if( cdsys_id>0 ) then
1082  if( iscoordneeds(g_localcoordsys(cdsys_id)) ) then
1083  coords=g_localcoordsys(cdsys_id)%CoordSys
1084  else
1085  ik=g_localcoordsys(cdsys_id)%node_ID(1)
1086  coords(1,:)= hecmesh%node(3*ik-2:3*ik)+fstrsolid%unode(3*ik-2:3*ik) &
1087  + fstrsolid%dunode(3*ik-2:3*ik)
1088  ik=g_localcoordsys(cdsys_id)%node_ID(2)
1089  coords(2,:)= hecmesh%node(3*ik-2:3*ik)+fstrsolid%unode(3*ik-2:3*ik) &
1090  + fstrsolid%dunode(3*ik-2:3*ik)
1091  ik=g_localcoordsys(cdsys_id)%node_ID(3)
1092  if(ik>0) coords(3,:)= hecmesh%node(3*ik-2:3*ik)+fstrsolid%unode(3*ik-2:3*ik) &
1093  + fstrsolid%dunode(3*ik-2:3*ik)
1094  endif
1095  endif
1096  end subroutine get_coordsys
1097 
1098  subroutine fstr_set_current_config_to_mesh(hecMESH,fstrSOLID,coord)
1099  implicit none
1100  type(hecmwst_local_mesh), intent(inout) :: hecMESH
1101  type (fstr_solid), intent(in) :: fstrSOLID
1102  real(kind=kreal), intent(out) :: coord(:)
1103  integer(kind=kint) :: i
1104  if(hecmesh%n_dof == 4) return
1105  do i = 1, hecmesh%nn_internal*min(hecmesh%n_dof,3)
1106  coord(i) = hecmesh%node(i)
1107  hecmesh%node(i) = coord(i)+fstrsolid%unode(i)+fstrsolid%dunode(i)
1108  enddo
1109  end subroutine fstr_set_current_config_to_mesh
1110 
1111  subroutine fstr_recover_initial_config_to_mesh(hecMESH,fstrSOLID,coord)
1112  implicit none
1113  type(hecmwst_local_mesh), intent(inout) :: hecMESH
1114  type (fstr_solid), intent(in) :: fstrSOLID
1115  real(kind=kreal), intent(in) :: coord(:)
1116  integer(kind=kint) :: i
1117  if(hecmesh%n_dof == 4) return
1118  do i = 1, hecmesh%nn_internal*min(hecmesh%n_dof,3)
1119  hecmesh%node(i) = coord(i)
1120  enddo
1122 
1123  subroutine fstr_solid_phys_zero(phys)
1124  implicit none
1125  type(fstr_solid_physic_val), pointer :: phys
1126  phys%STRAIN = 0.0d0
1127  phys%STRESS = 0.0d0
1128  phys%MISES = 0.0d0
1129  phys%ESTRAIN = 0.0d0
1130  phys%ESTRESS = 0.0d0
1131  phys%EMISES = 0.0d0
1132  phys%ENQM = 0.0d0
1133  end subroutine fstr_solid_phys_zero
1134 
1135  subroutine fstr_solid_phys_clear(fstrSOLID)
1136  implicit none
1137  type (fstr_solid) :: fstrSOLID
1138  integer(kind=kint) :: i
1139 
1140  if (associated(fstrsolid%SOLID)) then
1141  call fstr_solid_phys_zero(fstrsolid%SOLID)
1142  end if
1143  if (associated(fstrsolid%SHELL)) then
1144  call fstr_solid_phys_zero(fstrsolid%SHELL)
1145  do i=1,fstrsolid%max_lyr
1146  call fstr_solid_phys_zero(fstrsolid%SHELL%LAYER(i)%PLUS)
1147  call fstr_solid_phys_zero(fstrsolid%SHELL%LAYER(i)%MINUS)
1148  end do
1149  end if
1150  if (associated(fstrsolid%BEAM)) then
1151  call fstr_solid_phys_zero(fstrsolid%BEAM)
1152  end if
1153  end subroutine fstr_solid_phys_clear
1154 
1155  subroutine table_amp(hecMESH, fstrSOLID, cstep, jj_n_amp, time, f_t)
1156  type ( hecmwST_local_mesh ), intent(in) :: hecMESH
1157  type ( fstr_solid ), intent(in) :: fstrSOLID
1158  integer(kind=kint), intent(in) :: cstep
1159  integer(kind=kint), intent(in) :: jj_n_amp
1160  real(kind=kreal), intent(in) :: time
1161  real(kind=kreal), intent(out) :: f_t
1162 
1163  integer(kind=kint) :: i
1164  integer(kind=kint) :: jj1, jj2
1165  integer(kind=kint) :: s1, s2, flag
1166  real(kind=kreal) :: t_1, t_2, t_t, f_1, f_2, tincre
1167 
1168  s1 = 0; s2 = 0
1169 
1170  if( jj_n_amp <= 0 ) then ! Amplitude not defined
1171  if(myrank == 0) then
1172  write(imsg,*) 'internal error: amplitude table not defined'
1173  endif
1174  call hecmw_abort( hecmw_comm_get_comm() )
1175  endif
1176 
1177  tincre = fstrsolid%step_ctrl( cstep )%initdt
1178  jj1 = hecmesh%amp%amp_index(jj_n_amp - 1)
1179  jj2 = hecmesh%amp%amp_index(jj_n_amp)
1180 
1181  jj1 = jj1 + 2
1182  t_t = time-fstrsolid%step_ctrl(cstep)%starttime
1183 
1184  ! if(jj2 .eq. 0) then
1185  ! f_t = 1.0
1186  if(t_t .gt. hecmesh%amp%amp_table(jj2)) then
1187  f_t = hecmesh%amp%amp_val(jj2)
1188  else if(t_t .le. hecmesh%amp%amp_table(jj2)) then
1189  flag=0
1190  do i = jj1, jj2
1191  if(t_t .le. hecmesh%amp%amp_table(i)) then
1192  s2 = i
1193  s1 = i - 1
1194  flag = 1
1195  endif
1196  if( flag == 1 ) exit
1197  end do
1198 
1199  t_2 = hecmesh%amp%amp_table(s2)
1200  t_1 = hecmesh%amp%amp_table(s1)
1201  f_2 = hecmesh%amp%amp_val(s2)
1202  f_1 = hecmesh%amp%amp_val(s1)
1203  if( t_2-t_1 .lt. 1.0e-20) then
1204  if(myrank == 0) then
1205  write(imsg,*) 'stop due to t_2-t_1 <= 0'
1206  endif
1207  call hecmw_abort( hecmw_comm_get_comm())
1208  endif
1209  f_t = ((t_2*f_1 - t_1*f_2) + (f_2 - f_1)*t_t) / (t_2 - t_1)
1210  endif
1211 
1212  end subroutine table_amp
1213 
1214 end module m_fstr
m_fstr::ksmbicgstab
integer(kind=kint), parameter ksmbicgstab
Definition: m_fstr.f90:47
m_fstr::cntfilname
character(len=hecmw_filename_len) cntfilname
FILE NAME.
Definition: m_fstr.f90:103
m_fstr::kitrfloatingerror
integer(kind=kint), parameter kitrfloatingerror
Definition: m_fstr.f90:92
m_fstr::iecho
integer(kind=kint), pointer iecho
FLAG for ECHO/RESULT/POST.
Definition: m_fstr.f90:121
m_out
This module manages step information.
Definition: m_out.f90:6
m_fstr::restart_outlast
integer(kind=kint), parameter restart_outlast
restart type
Definition: m_fstr.f90:70
m_fstr::restartfilname
character(len=hecmw_filename_len) restartfilname
Definition: m_fstr.f90:104
m_step::iscontactactive
logical function iscontactactive(bnd, stepinfo)
Is contact condition in this step active.
Definition: m_step.f90:118
mcontactdef::tcontact
Structure to includes all info needed by contact calculation.
Definition: fstr_contact_def.F90:26
m_fstr::kon
integer(kind=kint), parameter kon
Definition: m_fstr.f90:32
m_fstr::kopss_material
integer(kind=kint), parameter kopss_material
Definition: m_fstr.f90:131
m_fstr::fstr_mat_init
subroutine fstr_mat_init(hecMAT)
Initializer of structure hecmwST_matrix.
Definition: m_fstr.f90:826
m_fstr::fstr_isloadactive
logical function fstr_isloadactive(fstrSOLID, nbc, cstep)
Definition: m_fstr.f90:1039
m_step::tparamautoinc
Definition: m_step.f90:51
m_fstr::kel361fi
integer(kind=kint), parameter kel361fi
Definition: m_fstr.f90:77
m_fstr::ksmgmres
integer(kind=kint), parameter ksmgmres
Definition: m_fstr.f90:48
m_fstr::ineutral
integer(kind=kint), pointer ineutral
Definition: m_fstr.f90:124
m_fstr::myrank
integer(kind=kint) myrank
PARALLEL EXECUTION.
Definition: m_fstr.f90:96
m_fstr::tsection
Data for section control.
Definition: m_fstr.f90:642
m_timepoint::time_points
Time points storage for output etc.
Definition: m_timepoint.f90:14
m_fstr::fstr_eigen
Package of data used by Lanczos eigenvalue solver.
Definition: m_fstr.f90:594
m_fstr::kststaticeigen
integer(kind=kint), parameter kststaticeigen
Definition: m_fstr.f90:42
m_fstr::itmax
integer(kind=kint) itmax
Definition: m_fstr.f90:141
m_common_struct::iscoordneeds
logical function iscoordneeds(coordsys)
if need to fetch global nodes' coordinate
Definition: m_common_struct.f90:62
mcontactdef
This module manages the data structure for contact calculation.
Definition: fstr_contact_def.F90:12
m_fstr::ksmcg
integer(kind=kint), parameter ksmcg
solver method (sm) !CAUTION : (<=100):indirect, (>100):direct
Definition: m_fstr.f90:46
m_fstr::fstr_freqanalysis_data
Definition: m_fstr.f90:585
m_fstr::fstr_heat
Data for HEAT ANSLYSIS (fstrHEAT)
Definition: m_fstr.f90:426
m_fstr::fstr_param_init
subroutine fstr_param_init(fstrPARAM, hecMESH)
Initializer of structure fstr_param.
Definition: m_fstr.f90:980
m_fstr::kel361bbar
integer(kind=kint), parameter kel361bbar
Definition: m_fstr.f90:78
m_fstr::ksmgpbicg
integer(kind=kint), parameter ksmgpbicg
Definition: m_fstr.f90:49
m_fstr::kbcinitial
integer(kind=kint), parameter kbcinitial
Definition: m_fstr.f90:66
m_fstr::kbctransit
integer(kind=kint), parameter kbctransit
Definition: m_fstr.f90:67
m_fstr::kbcfnastran
integer(kind=kint), parameter kbcfnastran
Definition: m_fstr.f90:64
m_fstr::kno
integer(kind=kint), parameter kno
Definition: m_fstr.f90:31
m_out::t_output_ctrl
output control such as output filename, output frequency etc.
Definition: m_out.f90:29
m_fstr::fstr_nullify_fstr_dynamic
subroutine fstr_nullify_fstr_dynamic(DY)
Definition: m_fstr.f90:783
m_fstr::fstr_solid
Definition: m_fstr.f90:239
m_common_struct::g_localcoordsys
type(tlocalcoordsys), dimension(:), pointer, save g_localcoordsys
Definition: m_common_struct.f90:19
m_fstr::fstr_dynamic
Data for DYNAMIC ANSLYSIS (fstrDYNAMIC)
Definition: m_fstr.f90:505
m_fstr::fstr_iscontactactive
logical function fstr_iscontactactive(fstrSOLID, nbc, cstep)
Definition: m_fstr.f90:1053
m_fstr::fstr_nullify_fstr_heat
subroutine fstr_nullify_fstr_heat(H)
Definition: m_fstr.f90:729
m_fstr::idbg
integer(kind=kint), parameter idbg
Definition: m_fstr.f90:111
m_fstr::fstr_param
FSTR INNER CONTROL PARAMETERS (fstrPARAM)
Definition: m_fstr.f90:154
m_fstr::fstr_nullify_fstr_eigen
subroutine fstr_nullify_fstr_eigen(E)
Definition: m_fstr.f90:808
m_fstr::kel341fi
integer(kind=kint), parameter kel341fi
section control
Definition: m_fstr.f90:74
m_fstr::opsstype
integer(kind=kint) opsstype
Definition: m_fstr.f90:132
m_fstr::knsmquasinewton
integer(kind=kint), parameter knsmquasinewton
Definition: m_fstr.f90:56
m_common_struct
This modules defines common structures for fem analysis.
Definition: m_common_struct.f90:6
m_fstr::kyes
integer(kind=kint), parameter kyes
CONSTANTS general.
Definition: m_fstr.f90:30
m_fstr::fstrpr
type(fstr_param), target fstrpr
GLOBAL VARIABLE INITIALIZED IN FSTR_SETUP.
Definition: m_fstr.f90:209
m_fstr::kitrcontinue
integer(kind=kint), parameter kitrcontinue
iteration control
Definition: m_fstr.f90:89
m_fstr::kitrconverged
integer(kind=kint), parameter kitrconverged
Definition: m_fstr.f90:90
m_fstr::tweldline
-1:not relation, >1:index of coupled_node
Definition: m_fstr.f90:629
m_fstr::eps
real(kind=kreal) eps
Definition: m_fstr.f90:142
mmechgauss
This modules defines a structure to record history dependent parameter in static analysis.
Definition: mechgauss.f90:6
m_fstr::kopss_solution
integer(kind=kint), parameter kopss_solution
Definition: m_fstr.f90:130
m_fstr::fstr_freqanalysis
Definition: m_fstr.f90:572
m_fstr::kel341sesns
integer(kind=kint), parameter kel341sesns
Definition: m_fstr.f90:75
m_fstr::kcaalagrange
integer(kind=kint), parameter kcaalagrange
Definition: m_fstr.f90:60
m_fstr::dt
real(kind=kreal) dt
ANALYSIS CONTROL for NLGEOM and HEAT.
Definition: m_fstr.f90:139
m_fstr::ifvs
integer(kind=kint), parameter ifvs
Definition: m_fstr.f90:112
m_fstr::fstr_isboundaryactive
logical function fstr_isboundaryactive(fstrSOLID, nbc, cstep)
Definition: m_fstr.f90:1029
fstr_sort_index
void fstr_sort_index(int *index_data, int *n)
Definition: fstr_sort_index.c:30
m_fstr::kcaslagrange
integer(kind=kint), parameter kcaslagrange
contact analysis algorithm
Definition: m_fstr.f90:59
m_fstr::kststatic
integer(kind=kint), parameter kststatic
Definition: m_fstr.f90:37
m_fstr::knsmnewton
integer(kind=kint), parameter knsmnewton
nonlinear solver method (nsm)
Definition: m_fstr.f90:55
m_fstr::iresout
integer(kind=kint), parameter iresout
Definition: m_fstr.f90:114
m_timepoint
This module manages timepoint information.
Definition: m_timepoint.f90:6
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
m_fstr::tinitialcondition
Definition: m_fstr.f90:144
m_fstr::nprocs
integer(kind=kint) nprocs
Definition: m_fstr.f90:97
m_fstr::ineu
integer(kind=kint), parameter ineu
Definition: m_fstr.f90:113
m_fstr::ivisual
integer(kind=kint), pointer ivisual
Definition: m_fstr.f90:123
m_fstr::ksmgmresr
integer(kind=kint), parameter ksmgmresr
Definition: m_fstr.f90:50
m_fstr::ref_temp
real(kind=kreal), pointer ref_temp
REFTEMP.
Definition: m_fstr.f90:136
m_fstr::g_initialcnd
type(tinitialcondition), dimension(:), pointer, save g_initialcnd
Definition: m_fstr.f90:151
m_fstr::fstr_couple
Data for coupling analysis.
Definition: m_fstr.f90:612
m_fstr::fstr_nullify_fstr_couple
subroutine fstr_nullify_fstr_couple(C)
Definition: m_fstr.f90:815
m_fstr::koff
integer(kind=kint), parameter koff
Definition: m_fstr.f90:33
m_fstr::kfloadtype_surf
integer(kind=kint), parameter kfloadtype_surf
Definition: m_fstr.f90:83
m_fstr::kfloadcase_re
integer(kind=kint), parameter kfloadcase_re
Definition: m_fstr.f90:85
m_fstr::fstr_solid_physic_val
Data for STATIC ANSLYSIS (fstrSOLID)
Definition: m_fstr.f90:212
m_fstr::fstr_nullify_fstr_freqanalysis
subroutine fstr_nullify_fstr_freqanalysis(f)
Definition: m_fstr.f90:795
m_fstr::kstheat
integer(kind=kint), parameter kstheat
Definition: m_fstr.f90:39
m_fstr::iresult
integer(kind=kint), pointer iresult
Definition: m_fstr.f90:122
m_fstr::restart_outall
integer(kind=kint), parameter restart_outall
Definition: m_fstr.f90:71
m_fstr::ksteigen
integer(kind=kint), parameter ksteigen
Definition: m_fstr.f90:38
m_fstr::ksmdirect
integer(kind=kint), parameter ksmdirect
Definition: m_fstr.f90:52
m_fstr::table_amp
subroutine table_amp(hecMESH, fstrSOLID, cstep, jj_n_amp, time, f_t)
Definition: m_fstr.f90:1156
m_fstr::fstr_nullify_fstr_param
subroutine fstr_nullify_fstr_param(P)
NULL POINTER SETTING TO AVOID RUNTIME ERROR.
Definition: m_fstr.f90:659
m_step::isboundaryactive
logical function isboundaryactive(bnd, stepinfo)
Is boundary condition in this step active.
Definition: m_step.f90:100
m_fstr::fstr_set_current_config_to_mesh
subroutine fstr_set_current_config_to_mesh(hecMESH, fstrSOLID, coord)
Definition: m_fstr.f90:1099
m_fstr::get_coordsys
subroutine get_coordsys(cdsys_ID, hecMESH, fstrSOLID, coords)
This subroutine fetch coords defined by local coordinate system.
Definition: m_fstr.f90:1074
hecmw
Definition: hecmw.f90:6
m_step::step_info
Step control such as active boundary condition, convergent condition etc.
Definition: m_step.f90:24
m_fstr::iwres
integer(kind=kint), pointer iwres
Definition: m_fstr.f90:126
m_fstr::ista
integer(kind=kint), parameter ista
Definition: m_fstr.f90:108
m_fstr::kstprecheck
integer(kind=kint), parameter kstprecheck
solution type (st)
Definition: m_fstr.f90:36
m_fstr::kfloadcase_im
integer(kind=kint), parameter kfloadcase_im
Definition: m_fstr.f90:86
m_fstr::kitrdiverged
integer(kind=kint), parameter kitrdiverged
Definition: m_fstr.f90:91
m_fstr::kel361fbar
integer(kind=kint), parameter kel361fbar
Definition: m_fstr.f90:80
m_fstr::kel361ic
integer(kind=kint), parameter kel361ic
Definition: m_fstr.f90:79
m_fstr::kfloadtype_node
integer(kind=kint), parameter kfloadtype_node
Definition: m_fstr.f90:82
m_fstr::kbcffstr
integer(kind=kint), parameter kbcffstr
boundary condition file type (bcf)
Definition: m_fstr.f90:63
m_fstr::hecmat_init
subroutine hecmat_init(hecMAT)
Definition: m_fstr.f90:866
m_fstr::iutb
integer(kind=kint), parameter iutb
Definition: m_fstr.f90:109
m_fstr::irres
integer(kind=kint), pointer irres
Definition: m_fstr.f90:125
m_fstr::nrres
integer(kind=kint), pointer nrres
Definition: m_fstr.f90:127
m_fstr::fstr_nullify_fstr_solid
subroutine fstr_nullify_fstr_solid(S)
Definition: m_fstr.f90:673
m_fstr::fstr_solid_phys_zero
subroutine fstr_solid_phys_zero(phys)
Definition: m_fstr.f90:1124
m_fstr::paracontactflag
logical paracontactflag
PARALLEL CONTACT FLAG.
Definition: m_fstr.f90:100
m_step::isloadactive
logical function isloadactive(bnd, stepinfo)
Is external load in this step active.
Definition: m_step.f90:109
m_fstr::sviarray
integer(kind=kint), dimension(100) sviarray
SOLVER CONTROL.
Definition: m_fstr.f90:117
m_fstr::etime
real(kind=kreal) etime
Definition: m_fstr.f90:140
m_fstr::ilog
integer(kind=kint), parameter ilog
FILE HANDLER.
Definition: m_fstr.f90:107
m_fstr::nprint
integer(kind=kint), pointer nprint
Definition: m_fstr.f90:128
m_fstr::svrarray
real(kind=kreal), dimension(100) svrarray
Definition: m_fstr.f90:118
mmechgauss::telement
All data should be recorded in every elements.
Definition: mechgauss.f90:32
m_fstr::hecmat_finalize
subroutine hecmat_finalize(hecMAT)
Definition: m_fstr.f90:923
m_step
This module manages step information.
Definition: m_step.f90:6
m_fstr::fstr_isembedactive
logical function fstr_isembedactive(fstrSOLID, nbc, cstep)
Definition: m_fstr.f90:1063
m_fstr::ksmgmresren
integer(kind=kint), parameter ksmgmresren
Definition: m_fstr.f90:51
m_fstr::fstr_recover_initial_config_to_mesh
subroutine fstr_recover_initial_config_to_mesh(hecMESH, fstrSOLID, coord)
Definition: m_fstr.f90:1112
m_fstr::fstr_solid_phys_clear
subroutine fstr_solid_phys_clear(fstrSOLID)
Definition: m_fstr.f90:1136
m_fstr::kstdynamic
integer(kind=kint), parameter kstdynamic
Definition: m_fstr.f90:40
m_fstr::kstnzprof
integer(kind=kint), parameter kstnzprof
Definition: m_fstr.f90:43
m_fstr::imsg
integer(kind=kint), parameter imsg
Definition: m_fstr.f90:110