FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
precheck_make_result.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 !-------------------------------------------------------------------------------
7 
8  implicit none
9 
10  private
11  public :: precheck_make_result
12  public :: precheck_write_result
13 
14 contains
15 
17  subroutine precheck_make_result(hecMESH, fstrRESULT, elem_vol, elem_asp)
18  use hecmw
19  use m_fstr
20  use elementinfo
21  implicit none
22  type(hecmwst_local_mesh), intent(in) :: hecmesh
23  type(hecmwst_result_data), intent(inout) :: fstrresult
24  real(kind=kreal), intent(in) :: elem_vol(:)
25  real(kind=kreal), intent(in) :: elem_asp(:)
26 
27  integer(kind=kint) :: n_ncomp, n_ecomp, n_nitem, n_eitem
28  integer(kind=kint) :: i, j, ie, ig, js, je, idx, isect, cid
29  integer(kind=kint) :: n_ngrp, n_egrp, n_sgrp
30  integer(kind=kint) :: ic, isurf, ic_type, stype, nn, j0, k
31  integer(kind=kint) :: snode(20)
32 
33  n_ngrp = hecmesh%node_group%n_grp
34  n_egrp = hecmesh%elem_group%n_grp
35  n_sgrp = hecmesh%surf_group%n_grp
36 
37  ! Node components: node groups + surf groups (as node flags)
38  n_ncomp = n_ngrp + n_sgrp
39  n_nitem = n_ngrp + n_sgrp ! each is scalar (dof=1)
40 
41  ! Element components: elem groups + surf groups + MATERIAL_ID + SECTION_ID + VOLUME + ASPECT_RATIO
42  n_ecomp = n_egrp + n_sgrp + 4
43  n_eitem = n_ecomp ! each is scalar (dof=1)
44 
45  ! Initialize result
46  fstrresult%ng_component = 0
47  fstrresult%nn_component = n_ncomp
48  fstrresult%ne_component = n_ecomp
49 
50  ! Allocate global (none)
51  allocate(fstrresult%ng_dof(0))
52  allocate(fstrresult%global_label(0))
53  allocate(fstrresult%global_val_item(0))
54 
55  ! Allocate node data
56  allocate(fstrresult%nn_dof(n_ncomp))
57  allocate(fstrresult%node_label(n_ncomp))
58  allocate(fstrresult%node_val_item(n_nitem * hecmesh%n_node))
59  fstrresult%nn_dof = 1
60  fstrresult%node_val_item = 0.0d0
61 
62  ! Allocate element data
63  allocate(fstrresult%ne_dof(n_ecomp))
64  allocate(fstrresult%elem_label(n_ecomp))
65  allocate(fstrresult%elem_val_item(n_eitem * hecmesh%n_elem))
66  fstrresult%ne_dof = 1
67  fstrresult%elem_val_item = 0.0d0
68 
69  ! --- Fill node group flags ---
70  do ig = 1, n_ngrp
71  fstrresult%node_label(ig) = 'NGRP_' // trim(hecmesh%node_group%grp_name(ig))
72  js = hecmesh%node_group%grp_index(ig-1)
73  je = hecmesh%node_group%grp_index(ig)
74  do j = js+1, je
75  idx = hecmesh%node_group%grp_item(j)
76  if(idx >= 1 .and. idx <= hecmesh%n_node) then
77  fstrresult%node_val_item(n_nitem*(idx-1) + ig) = 1.0d0
78  endif
79  enddo
80  enddo
81 
82  ! --- Fill surface group node flags ---
83  do ig = 1, n_sgrp
84  fstrresult%node_label(n_ngrp + ig) = 'SGRP_' // trim(hecmesh%surf_group%grp_name(ig))
85  js = hecmesh%surf_group%grp_index(ig-1)
86  je = hecmesh%surf_group%grp_index(ig)
87  do j = js+1, je
88  ic = hecmesh%surf_group%grp_item(2*j-1) ! element ID
89  isurf = hecmesh%surf_group%grp_item(2*j) ! face number
90  if(ic < 1 .or. ic > hecmesh%n_elem) cycle
91  ic_type = hecmesh%elem_type(ic)
92  call getsubface(ic_type, isurf, stype, snode)
93  nn = getnumberofnodes(stype)
94  j0 = hecmesh%elem_node_index(ic-1)
95  do k = 1, nn
96  idx = hecmesh%elem_node_item(j0 + snode(k))
97  if(idx >= 1 .and. idx <= hecmesh%n_node) then
98  fstrresult%node_val_item(n_nitem*(idx-1) + n_ngrp + ig) = 1.0d0
99  endif
100  enddo
101  enddo
102  enddo
103 
104  ! --- Fill element group flags ---
105  do ig = 1, n_egrp
106  fstrresult%elem_label(ig) = 'EGRP_' // trim(hecmesh%elem_group%grp_name(ig))
107  js = hecmesh%elem_group%grp_index(ig-1)
108  je = hecmesh%elem_group%grp_index(ig)
109  do j = js+1, je
110  idx = hecmesh%elem_group%grp_item(j)
111  if(idx >= 1 .and. idx <= hecmesh%n_elem) then
112  fstrresult%elem_val_item(n_eitem*(idx-1) + ig) = 1.0d0
113  endif
114  enddo
115  enddo
116 
117  ! --- Fill surface group (element with face number) ---
118  do ig = 1, n_sgrp
119  fstrresult%elem_label(n_egrp + ig) = 'SGRP_FACE_' // trim(hecmesh%surf_group%grp_name(ig))
120  js = hecmesh%surf_group%grp_index(ig-1)
121  je = hecmesh%surf_group%grp_index(ig)
122  do j = js+1, je
123  idx = hecmesh%surf_group%grp_item(2*j-1) ! element ID
124  if(idx >= 1 .and. idx <= hecmesh%n_elem) then
125  fstrresult%elem_val_item(n_eitem*(idx-1) + n_egrp + ig) = &
126  dble(hecmesh%surf_group%grp_item(2*j)) ! face number
127  endif
128  enddo
129  enddo
130 
131  ! --- Fill MATERIAL_ID ---
132  i = n_egrp + n_sgrp + 1
133  fstrresult%elem_label(i) = 'MATERIAL_ID'
134  do ie = 1, hecmesh%n_elem
135  isect = hecmesh%section_ID(ie)
136  if(isect >= 1 .and. isect <= hecmesh%section%n_sect) then
137  cid = hecmesh%section%sect_mat_ID_item(isect)
138  fstrresult%elem_val_item(n_eitem*(ie-1) + i) = dble(cid)
139  endif
140  enddo
141 
142  ! --- Fill SECTION_ID ---
143  i = n_egrp + n_sgrp + 2
144  fstrresult%elem_label(i) = 'SECTION_ID'
145  do ie = 1, hecmesh%n_elem
146  fstrresult%elem_val_item(n_eitem*(ie-1) + i) = dble(hecmesh%section_ID(ie))
147  enddo
148 
149  ! --- Fill VOLUME ---
150  i = n_egrp + n_sgrp + 3
151  fstrresult%elem_label(i) = 'VOLUME'
152  do ie = 1, hecmesh%n_elem
153  fstrresult%elem_val_item(n_eitem*(ie-1) + i) = elem_vol(ie)
154  enddo
155 
156  ! --- Fill ASPECT_RATIO ---
157  i = n_egrp + n_sgrp + 4
158  fstrresult%elem_label(i) = 'ASPECT_RATIO'
159  do ie = 1, hecmesh%n_elem
160  fstrresult%elem_val_item(n_eitem*(ie-1) + i) = elem_asp(ie)
161  enddo
162 
163  end subroutine precheck_make_result
164 
166  subroutine precheck_write_result(hecMESH, elem_vol, elem_asp)
167  use hecmw
168  use m_fstr
169  use hecmw_result
170  use elementinfo
171  implicit none
172  type(hecmwst_local_mesh), intent(in) :: hecmesh
173  real(kind=kreal), intent(in) :: elem_vol(:)
174  real(kind=kreal), intent(in) :: elem_asp(:)
175 
176  integer(kind=kint) :: ig, j, js, je, idx, ie, isect, cid
177  integer(kind=kint) :: n_ngrp, n_egrp, n_sgrp
178  integer(kind=kint) :: ic, isurf, ic_type, stype, nn, j0, k
179  integer(kind=kint) :: snode(20)
180  character(len=HECMW_HEADER_LEN) :: header
181  character(len=HECMW_MSG_LEN) :: comment
182  character(len=HECMW_NAME_LEN) :: label, nameid, addfname
183  real(kind=kreal), allocatable :: work_n(:), work_e(:)
184 
185  allocate(work_n(hecmesh%n_node))
186  allocate(work_e(hecmesh%n_elem))
187 
188  header = '*fstrresult'
189  comment = 'elemcheck_result'
190  call hecmw_result_init(hecmesh, 0, header, comment)
191 
192  n_ngrp = hecmesh%node_group%n_grp
193  n_egrp = hecmesh%elem_group%n_grp
194  n_sgrp = hecmesh%surf_group%n_grp
195 
196  ! --- Node groups ---
197  do ig = 1, n_ngrp
198  work_n = 0.0d0
199  js = hecmesh%node_group%grp_index(ig-1)
200  je = hecmesh%node_group%grp_index(ig)
201  do j = js+1, je
202  idx = hecmesh%node_group%grp_item(j)
203  if(idx >= 1 .and. idx <= hecmesh%n_node) work_n(idx) = 1.0d0
204  enddo
205  label = 'NGRP_' // trim(hecmesh%node_group%grp_name(ig))
206  call hecmw_result_add(hecmw_result_dtype_node, 1, label, work_n)
207  enddo
208 
209  ! --- Surface groups (node flags) ---
210  do ig = 1, n_sgrp
211  work_n = 0.0d0
212  js = hecmesh%surf_group%grp_index(ig-1)
213  je = hecmesh%surf_group%grp_index(ig)
214  do j = js+1, je
215  ic = hecmesh%surf_group%grp_item(2*j-1)
216  isurf = hecmesh%surf_group%grp_item(2*j)
217  if(ic < 1 .or. ic > hecmesh%n_elem) cycle
218  ic_type = hecmesh%elem_type(ic)
219  call getsubface(ic_type, isurf, stype, snode)
220  nn = getnumberofnodes(stype)
221  j0 = hecmesh%elem_node_index(ic-1)
222  do k = 1, nn
223  idx = hecmesh%elem_node_item(j0 + snode(k))
224  if(idx >= 1 .and. idx <= hecmesh%n_node) work_n(idx) = 1.0d0
225  enddo
226  enddo
227  label = 'SGRP_' // trim(hecmesh%surf_group%grp_name(ig))
228  call hecmw_result_add(hecmw_result_dtype_node, 1, label, work_n)
229  enddo
230 
231  ! --- Element groups ---
232  do ig = 1, n_egrp
233  work_e = 0.0d0
234  js = hecmesh%elem_group%grp_index(ig-1)
235  je = hecmesh%elem_group%grp_index(ig)
236  do j = js+1, je
237  idx = hecmesh%elem_group%grp_item(j)
238  if(idx >= 1 .and. idx <= hecmesh%n_elem) work_e(idx) = 1.0d0
239  enddo
240  label = 'EGRP_' // trim(hecmesh%elem_group%grp_name(ig))
241  call hecmw_result_add(hecmw_result_dtype_elem, 1, label, work_e)
242  enddo
243 
244  ! --- Surface groups ---
245  do ig = 1, n_sgrp
246  work_e = 0.0d0
247  js = hecmesh%surf_group%grp_index(ig-1)
248  je = hecmesh%surf_group%grp_index(ig)
249  do j = js+1, je
250  idx = hecmesh%surf_group%grp_item(2*j-1)
251  if(idx >= 1 .and. idx <= hecmesh%n_elem) then
252  work_e(idx) = dble(hecmesh%surf_group%grp_item(2*j))
253  endif
254  enddo
255  label = 'SGRP_FACE_' // trim(hecmesh%surf_group%grp_name(ig))
256  call hecmw_result_add(hecmw_result_dtype_elem, 1, label, work_e)
257  enddo
258 
259  ! --- MATERIAL_ID ---
260  do ie = 1, hecmesh%n_elem
261  isect = hecmesh%section_ID(ie)
262  if(isect >= 1 .and. isect <= hecmesh%section%n_sect) then
263  cid = hecmesh%section%sect_mat_ID_item(isect)
264  work_e(ie) = dble(cid)
265  else
266  work_e(ie) = 0.0d0
267  endif
268  enddo
269  call hecmw_result_add(hecmw_result_dtype_elem, 1, 'MATERIAL_ID', work_e)
270 
271  ! --- SECTION_ID ---
272  do ie = 1, hecmesh%n_elem
273  work_e(ie) = dble(hecmesh%section_ID(ie))
274  enddo
275  call hecmw_result_add(hecmw_result_dtype_elem, 1, 'SECTION_ID', work_e)
276 
277  ! --- VOLUME ---
278  call hecmw_result_add(hecmw_result_dtype_elem, 1, 'VOLUME', elem_vol)
279 
280  ! --- ASPECT_RATIO ---
281  call hecmw_result_add(hecmw_result_dtype_elem, 1, 'ASPECT_RATIO', elem_asp)
282 
283  ! --- WRITE ---
284  nameid = 'fstrRES'
285  addfname = '_precheck'
286  call hecmw_result_write_by_addfname(nameid, addfname)
288 
289  deallocate(work_n)
290  deallocate(work_e)
291 
292  end subroutine precheck_write_result
293 
294 end module m_precheck_make_result
This module encapsulate the basic functions of all elements provide by this software.
Definition: element.f90:43
integer(kind=kind(2)) function getnumberofnodes(etype)
Obtain number of nodes of the element.
Definition: element.f90:131
subroutine getsubface(intype, innumber, outtype, nodes)
Find the definition of surface of the element.
Definition: element.f90:193
I/O and Utility.
subroutine, public hecmw_result_write_by_addfname(name_ID, addfname)
subroutine, public hecmw_result_add(dtype, n_dof, label, data)
subroutine, public hecmw_result_finalize()
subroutine, public hecmw_result_init(hecMESH, i_step, header, comment)
integer(kind=kint), parameter, public hecmw_result_dtype_node
integer(kind=kint), parameter, public hecmw_result_dtype_elem
Definition: hecmw.f90:6
This module defines common data and basic structures for analysis.
Definition: m_fstr.F90:15
This module constructs result data for ELEMCHECK visualization/result output.
subroutine, public precheck_make_result(hecMESH, fstrRESULT, elem_vol, elem_asp)
Build hecmwST_result_data for visualization via hecmw_visualize.
subroutine, public precheck_write_result(hecMESH, elem_vol, elem_asp)
Write check results to result file via hecmw_result_add/write interface.