FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hecmw_solver_scaling_66.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 
7  use hecmw_util
10  implicit none
11 
12  private
13  real(kind=kreal), private, allocatable :: scale(:)
14 
17 
18 contains
19 
20  subroutine hecmw_solver_scaling_fw_66(hecMESH, hecMAT, COMMtime)
21  implicit none
22  type (hecmwst_local_mesh), intent(in) :: hecmesh
23  type (hecmwst_matrix), intent(inout) :: hecmat
24  real(kind=kreal), intent(inout) :: commtime
25  integer(kind=kint) :: n, np, ndof
26  real(kind=kreal), pointer :: d(:), al(:), au(:), b(:)
27  integer(kind=kint), pointer :: inl(:), ial(:), inu(:), iau(:)
28  integer(kind=kint) :: i, k, ip1, ip2, ip3, ip4, ip5, ip6
29  integer(kind=kint) :: iq1, iq2, iq3, iq4, iq5, iq6
30  integer(kind=kint) :: isl, iel, isu, ieu, inod
31  real(kind=kreal) :: start_time, end_time
32 
33  if (hecmw_mat_get_scaling(hecmat).eq.0) return
34 
35  n = hecmat%N
36  np = hecmat%NP
37  ndof = hecmat%NDOF
38  d => hecmat%D
39  al => hecmat%AL
40  au => hecmat%AU
41  inl => hecmat%indexL
42  ial => hecmat%itemL
43  inu => hecmat%indexU
44  iau => hecmat%itemU
45  b => hecmat%B
46 
47  allocate(scale(ndof*np))
48 
49  do i= 1, n
50  scale(6*i-5)= 1.d0/dsqrt(dabs(d(36*i-35)))
51  scale(6*i-4)= 1.d0/dsqrt(dabs(d(36*i-28)))
52  scale(6*i-3)= 1.d0/dsqrt(dabs(d(36*i-21)))
53  scale(6*i-2)= 1.d0/dsqrt(dabs(d(36*i-14)))
54  scale(6*i-1)= 1.d0/dsqrt(dabs(d(36*i-7)))
55  scale(6*i )= 1.d0/dsqrt(dabs(d(36*i )))
56  enddo
57 
58  start_time= hecmw_wtime()
59  call hecmw_update_r (hecmesh, scale, hecmesh%n_node, ndof)
60  end_time= hecmw_wtime()
61  commtime = commtime + end_time - start_time
62 
63  do i= 1, np
64  ip1= 6*i-5
65  ip2= 6*i-4
66  ip3= 6*i-3
67  ip4= 6*i-2
68  ip5= 6*i-1
69  ip6= 6*i
70 
71  d(36*i-35)= d(36*i-35)*scale(ip1)*scale(ip1)
72  d(36*i-34)= d(36*i-34)*scale(ip1)*scale(ip2)
73  d(36*i-33)= d(36*i-33)*scale(ip1)*scale(ip3)
74  d(36*i-32)= d(36*i-32)*scale(ip1)*scale(ip4)
75  d(36*i-31)= d(36*i-31)*scale(ip1)*scale(ip5)
76  d(36*i-30)= d(36*i-30)*scale(ip1)*scale(ip6)
77 
78  d(36*i-29)= d(36*i-29)*scale(ip2)*scale(ip1)
79  d(36*i-28)= d(36*i-28)*scale(ip2)*scale(ip2)
80  d(36*i-27)= d(36*i-27)*scale(ip2)*scale(ip3)
81  d(36*i-26)= d(36*i-26)*scale(ip2)*scale(ip4)
82  d(36*i-25)= d(36*i-25)*scale(ip2)*scale(ip5)
83  d(36*i-24)= d(36*i-24)*scale(ip2)*scale(ip6)
84 
85  d(36*i-23)= d(36*i-23)*scale(ip3)*scale(ip1)
86  d(36*i-22)= d(36*i-22)*scale(ip3)*scale(ip2)
87  d(36*i-21)= d(36*i-21)*scale(ip3)*scale(ip3)
88  d(36*i-20)= d(36*i-20)*scale(ip3)*scale(ip4)
89  d(36*i-19)= d(36*i-19)*scale(ip3)*scale(ip5)
90  d(36*i-18)= d(36*i-18)*scale(ip3)*scale(ip6)
91 
92  d(36*i-17)= d(36*i-17)*scale(ip4)*scale(ip1)
93  d(36*i-16)= d(36*i-16)*scale(ip4)*scale(ip2)
94  d(36*i-15)= d(36*i-15)*scale(ip4)*scale(ip3)
95  d(36*i-14)= d(36*i-14)*scale(ip4)*scale(ip4)
96  d(36*i-13)= d(36*i-13)*scale(ip4)*scale(ip5)
97  d(36*i-12)= d(36*i-12)*scale(ip4)*scale(ip6)
98 
99  d(36*i-11)= d(36*i-11)*scale(ip5)*scale(ip1)
100  d(36*i-10)= d(36*i-10)*scale(ip5)*scale(ip2)
101  d(36*i-9 )= d(36*i-9 )*scale(ip5)*scale(ip3)
102  d(36*i-8 )= d(36*i-8 )*scale(ip5)*scale(ip4)
103  d(36*i-7 )= d(36*i-7 )*scale(ip5)*scale(ip5)
104  d(36*i-6 )= d(36*i-6 )*scale(ip5)*scale(ip6)
105 
106  d(36*i-5 )= d(36*i-5 )*scale(ip6)*scale(ip1)
107  d(36*i-4 )= d(36*i-4 )*scale(ip6)*scale(ip2)
108  d(36*i-3 )= d(36*i-3 )*scale(ip6)*scale(ip3)
109  d(36*i-2 )= d(36*i-2 )*scale(ip6)*scale(ip4)
110  d(36*i-1 )= d(36*i-1 )*scale(ip6)*scale(ip5)
111  d(36*i )= d(36*i )*scale(ip6)*scale(ip6)
112 
113  isl= inl(i-1) + 1
114  iel= inl(i )
115  !*voption indep (IAL,AL,SCALE)
116  do k= isl, iel
117  inod= ial(k)
118  iq1= 6*inod - 5
119  iq2= 6*inod - 4
120  iq3= 6*inod - 3
121  iq4= 6*inod - 2
122  iq5= 6*inod - 1
123  iq6= 6*inod
124 
125  al(36*k-35)= al(36*k-35)*scale(ip1)*scale(iq1)
126  al(36*k-34)= al(36*k-34)*scale(ip1)*scale(iq2)
127  al(36*k-33)= al(36*k-33)*scale(ip1)*scale(iq3)
128  al(36*k-32)= al(36*k-32)*scale(ip1)*scale(iq4)
129  al(36*k-31)= al(36*k-31)*scale(ip1)*scale(iq5)
130  al(36*k-30)= al(36*k-30)*scale(ip1)*scale(iq6)
131 
132  al(36*k-29)= al(36*k-29)*scale(ip2)*scale(iq1)
133  al(36*k-28)= al(36*k-28)*scale(ip2)*scale(iq2)
134  al(36*k-27)= al(36*k-27)*scale(ip2)*scale(iq3)
135  al(36*k-26)= al(36*k-26)*scale(ip2)*scale(iq4)
136  al(36*k-25)= al(36*k-25)*scale(ip2)*scale(iq5)
137  al(36*k-24)= al(36*k-24)*scale(ip2)*scale(iq6)
138 
139  al(36*k-23)= al(36*k-23)*scale(ip3)*scale(iq1)
140  al(36*k-22)= al(36*k-22)*scale(ip3)*scale(iq2)
141  al(36*k-21)= al(36*k-21)*scale(ip3)*scale(iq3)
142  al(36*k-20)= al(36*k-20)*scale(ip3)*scale(iq4)
143  al(36*k-19)= al(36*k-19)*scale(ip3)*scale(iq5)
144  al(36*k-18)= al(36*k-18)*scale(ip3)*scale(iq6)
145 
146  al(36*k-17)= al(36*k-17)*scale(ip4)*scale(iq1)
147  al(36*k-16)= al(36*k-16)*scale(ip4)*scale(iq2)
148  al(36*k-15)= al(36*k-15)*scale(ip4)*scale(iq3)
149  al(36*k-14)= al(36*k-14)*scale(ip4)*scale(iq4)
150  al(36*k-13)= al(36*k-13)*scale(ip4)*scale(iq5)
151  al(36*k-12)= al(36*k-12)*scale(ip4)*scale(iq6)
152 
153  al(36*k-11)= al(36*k-11)*scale(ip5)*scale(iq1)
154  al(36*k-10)= al(36*k-10)*scale(ip5)*scale(iq2)
155  al(36*k-9 )= al(36*k-9 )*scale(ip5)*scale(iq3)
156  al(36*k-8 )= al(36*k-8 )*scale(ip5)*scale(iq4)
157  al(36*k-7 )= al(36*k-7 )*scale(ip5)*scale(iq5)
158  al(36*k-6 )= al(36*k-6 )*scale(ip5)*scale(iq6)
159 
160  al(36*k-5 )= al(36*k-5 )*scale(ip6)*scale(iq1)
161  al(36*k-4 )= al(36*k-4 )*scale(ip6)*scale(iq2)
162  al(36*k-3 )= al(36*k-3 )*scale(ip6)*scale(iq3)
163  al(36*k-2 )= al(36*k-2 )*scale(ip6)*scale(iq4)
164  al(36*k-1 )= al(36*k-1 )*scale(ip6)*scale(iq5)
165  al(36*k )= al(36*k )*scale(ip6)*scale(iq6)
166  enddo
167 
168  isu= inu(i-1) + 1
169  ieu= inu(i )
170  !*voption indep (IAU,AU,SCALE)
171  do k= isu, ieu
172  inod= iau(k)
173  iq1= 6*inod - 5
174  iq2= 6*inod - 4
175  iq3= 6*inod - 3
176  iq4= 6*inod - 2
177  iq5= 6*inod - 1
178  iq6= 6*inod
179 
180  au(36*k-35)= au(36*k-35)*scale(ip1)*scale(iq1)
181  au(36*k-34)= au(36*k-34)*scale(ip1)*scale(iq2)
182  au(36*k-33)= au(36*k-33)*scale(ip1)*scale(iq3)
183  au(36*k-32)= au(36*k-32)*scale(ip1)*scale(iq4)
184  au(36*k-31)= au(36*k-31)*scale(ip1)*scale(iq5)
185  au(36*k-30)= au(36*k-30)*scale(ip1)*scale(iq6)
186 
187  au(36*k-29)= au(36*k-29)*scale(ip2)*scale(iq1)
188  au(36*k-28)= au(36*k-28)*scale(ip2)*scale(iq2)
189  au(36*k-27)= au(36*k-27)*scale(ip2)*scale(iq3)
190  au(36*k-26)= au(36*k-26)*scale(ip2)*scale(iq4)
191  au(36*k-25)= au(36*k-25)*scale(ip2)*scale(iq5)
192  au(36*k-24)= au(36*k-24)*scale(ip2)*scale(iq6)
193 
194  au(36*k-23)= au(36*k-23)*scale(ip3)*scale(iq1)
195  au(36*k-22)= au(36*k-22)*scale(ip3)*scale(iq2)
196  au(36*k-21)= au(36*k-21)*scale(ip3)*scale(iq3)
197  au(36*k-20)= au(36*k-20)*scale(ip3)*scale(iq4)
198  au(36*k-19)= au(36*k-19)*scale(ip3)*scale(iq5)
199  au(36*k-18)= au(36*k-18)*scale(ip3)*scale(iq6)
200 
201  au(36*k-17)= au(36*k-17)*scale(ip4)*scale(iq1)
202  au(36*k-16)= au(36*k-16)*scale(ip4)*scale(iq2)
203  au(36*k-15)= au(36*k-15)*scale(ip4)*scale(iq3)
204  au(36*k-14)= au(36*k-14)*scale(ip4)*scale(iq4)
205  au(36*k-13)= au(36*k-13)*scale(ip4)*scale(iq5)
206  au(36*k-12)= au(36*k-12)*scale(ip4)*scale(iq6)
207 
208  au(36*k-11)= au(36*k-11)*scale(ip5)*scale(iq1)
209  au(36*k-10)= au(36*k-10)*scale(ip5)*scale(iq2)
210  au(36*k-9 )= au(36*k-9 )*scale(ip5)*scale(iq3)
211  au(36*k-8 )= au(36*k-8 )*scale(ip5)*scale(iq4)
212  au(36*k-7 )= au(36*k-7 )*scale(ip5)*scale(iq5)
213  au(36*k-6 )= au(36*k-6 )*scale(ip5)*scale(iq6)
214 
215  au(36*k-5 )= au(36*k-5 )*scale(ip6)*scale(iq1)
216  au(36*k-4 )= au(36*k-4 )*scale(ip6)*scale(iq2)
217  au(36*k-3 )= au(36*k-3 )*scale(ip6)*scale(iq3)
218  au(36*k-2 )= au(36*k-2 )*scale(ip6)*scale(iq4)
219  au(36*k-1 )= au(36*k-1 )*scale(ip6)*scale(iq5)
220  au(36*k )= au(36*k )*scale(ip6)*scale(iq6)
221  enddo
222  enddo
223  !*voption indep (B,SCALE)
224  do i= 1, n
225  b(6*i-5)= b(6*i-5) * scale(6*i-5)
226  b(6*i-4)= b(6*i-4) * scale(6*i-4)
227  b(6*i-3)= b(6*i-3) * scale(6*i-3)
228  b(6*i-2)= b(6*i-2) * scale(6*i-2)
229  b(6*i-1)= b(6*i-1) * scale(6*i-1)
230  b(6*i )= b(6*i ) * scale(6*i )
231  enddo
232  end subroutine hecmw_solver_scaling_fw_66
233 
234  subroutine hecmw_solver_scaling_bk_66(hecMAT)
236  implicit none
237  type (hecmwst_matrix), intent(inout) :: hecmat
238  integer(kind=kint) :: n, np, ndof
239  real(kind=kreal), pointer :: d(:), al(:), au(:), b(:), x(:)
240  integer(kind=kint), pointer :: inl(:), ial(:), inu(:), iau(:)
241  integer(kind=kint) :: i, k, ip1, ip2, ip3, ip4, ip5, ip6
242  integer(kind=kint) :: iq1, iq2, iq3, iq4, iq5, iq6
243  integer(kind=kint) :: isl, iel, isu, ieu, inod
244 
245  if (hecmw_mat_get_scaling(hecmat).eq.0) return
246 
247  n = hecmat%N
248  np = hecmat%NP
249  ndof = hecmat%NDOF
250  d => hecmat%D
251  al => hecmat%AL
252  au => hecmat%AU
253  inl => hecmat%indexL
254  ial => hecmat%itemL
255  inu => hecmat%indexU
256  iau => hecmat%itemU
257  b => hecmat%B
258  x => hecmat%X
259 
260  !*voption indep (X,B,SCALE)
261  do i= 1, n
262  x(6*i-5)= x(6*i-5) * scale(6*i-5)
263  x(6*i-4)= x(6*i-4) * scale(6*i-4)
264  x(6*i-3)= x(6*i-3) * scale(6*i-3)
265  x(6*i-2)= x(6*i-2) * scale(6*i-2)
266  x(6*i-1)= x(6*i-1) * scale(6*i-1)
267  x(6*i )= x(6*i ) * scale(6*i )
268 
269  b(6*i-5)= b(6*i-5) / scale(6*i-5)
270  b(6*i-4)= b(6*i-4) / scale(6*i-4)
271  b(6*i-3)= b(6*i-3) / scale(6*i-3)
272  b(6*i-2)= b(6*i-2) / scale(6*i-2)
273  b(6*i-1)= b(6*i-1) / scale(6*i-1)
274  b(6*i )= b(6*i ) / scale(6*i )
275  enddo
276 
277  do i= 1, np
278  ip1= 6*i-5
279  ip2= 6*i-4
280  ip3= 6*i-3
281  ip4= 6*i-2
282  ip5= 6*i-1
283  ip6= 6*i
284 
285  d(36*i-35)= d(36*i-35)/(scale(ip1)*scale(ip1))
286  d(36*i-34)= d(36*i-34)/(scale(ip1)*scale(ip2))
287  d(36*i-33)= d(36*i-33)/(scale(ip1)*scale(ip3))
288  d(36*i-32)= d(36*i-32)/(scale(ip1)*scale(ip4))
289  d(36*i-31)= d(36*i-31)/(scale(ip1)*scale(ip5))
290  d(36*i-30)= d(36*i-30)/(scale(ip1)*scale(ip6))
291 
292  d(36*i-29)= d(36*i-29)/(scale(ip2)*scale(ip1))
293  d(36*i-28)= d(36*i-28)/(scale(ip2)*scale(ip2))
294  d(36*i-27)= d(36*i-27)/(scale(ip2)*scale(ip3))
295  d(36*i-26)= d(36*i-26)/(scale(ip2)*scale(ip4))
296  d(36*i-25)= d(36*i-25)/(scale(ip2)*scale(ip5))
297  d(36*i-24)= d(36*i-24)/(scale(ip2)*scale(ip6))
298 
299  d(36*i-23)= d(36*i-23)/(scale(ip3)*scale(ip1))
300  d(36*i-22)= d(36*i-22)/(scale(ip3)*scale(ip2))
301  d(36*i-21)= d(36*i-21)/(scale(ip3)*scale(ip3))
302  d(36*i-20)= d(36*i-20)/(scale(ip3)*scale(ip4))
303  d(36*i-19)= d(36*i-19)/(scale(ip3)*scale(ip5))
304  d(36*i-18)= d(36*i-18)/(scale(ip3)*scale(ip6))
305 
306  d(36*i-17)= d(36*i-17)/(scale(ip4)*scale(ip1))
307  d(36*i-16)= d(36*i-16)/(scale(ip4)*scale(ip2))
308  d(36*i-15)= d(36*i-15)/(scale(ip4)*scale(ip3))
309  d(36*i-14)= d(36*i-14)/(scale(ip4)*scale(ip4))
310  d(36*i-13)= d(36*i-13)/(scale(ip4)*scale(ip5))
311  d(36*i-12)= d(36*i-12)/(scale(ip4)*scale(ip6))
312 
313  d(36*i-11)= d(36*i-11)/(scale(ip5)*scale(ip1))
314  d(36*i-10)= d(36*i-10)/(scale(ip5)*scale(ip2))
315  d(36*i-9 )= d(36*i-9 )/(scale(ip5)*scale(ip3))
316  d(36*i-8 )= d(36*i-8 )/(scale(ip5)*scale(ip4))
317  d(36*i-7 )= d(36*i-7 )/(scale(ip5)*scale(ip5))
318  d(36*i-6 )= d(36*i-6 )/(scale(ip5)*scale(ip6))
319 
320  d(36*i-5 )= d(36*i-5 )/(scale(ip6)*scale(ip1))
321  d(36*i-4 )= d(36*i-4 )/(scale(ip6)*scale(ip2))
322  d(36*i-3 )= d(36*i-3 )/(scale(ip6)*scale(ip3))
323  d(36*i-2 )= d(36*i-2 )/(scale(ip6)*scale(ip4))
324  d(36*i-1 )= d(36*i-1 )/(scale(ip6)*scale(ip5))
325  d(36*i )= d(36*i )/(scale(ip6)*scale(ip6))
326 
327  isl= inl(i-1) + 1
328  iel= inl(i )
329  !*voption indep (IAL,AL,SCALE)
330  do k= isl, iel
331  inod= ial(k)
332  iq1= 6*inod - 5
333  iq2= 6*inod - 5
334  iq3= 6*inod - 3
335  iq4= 6*inod - 2
336  iq5= 6*inod - 1
337  iq6= 6*inod
338 
339  al(36*k-35)= al(36*k-35)/(scale(ip1)*scale(iq1))
340  al(36*k-34)= al(36*k-34)/(scale(ip1)*scale(iq2))
341  al(36*k-33)= al(36*k-33)/(scale(ip1)*scale(iq3))
342  al(36*k-32)= al(36*k-32)/(scale(ip1)*scale(iq4))
343  al(36*k-31)= al(36*k-31)/(scale(ip1)*scale(iq5))
344  al(36*k-30)= al(36*k-30)/(scale(ip1)*scale(iq6))
345 
346  al(36*k-29)= al(36*k-29)/(scale(ip2)*scale(iq1))
347  al(36*k-28)= al(36*k-28)/(scale(ip2)*scale(iq2))
348  al(36*k-27)= al(36*k-27)/(scale(ip2)*scale(iq3))
349  al(36*k-26)= al(36*k-26)/(scale(ip2)*scale(iq4))
350  al(36*k-25)= al(36*k-25)/(scale(ip2)*scale(iq5))
351  al(36*k-24)= al(36*k-24)/(scale(ip2)*scale(iq6))
352 
353  al(36*k-23)= al(36*k-23)/(scale(ip3)*scale(iq1))
354  al(36*k-22)= al(36*k-22)/(scale(ip3)*scale(iq2))
355  al(36*k-21)= al(36*k-21)/(scale(ip3)*scale(iq3))
356  al(36*k-20)= al(36*k-20)/(scale(ip3)*scale(iq4))
357  al(36*k-19)= al(36*k-19)/(scale(ip3)*scale(iq5))
358  al(36*k-18)= al(36*k-18)/(scale(ip3)*scale(iq6))
359 
360  al(36*k-17)= al(36*k-17)/(scale(ip4)*scale(iq1))
361  al(36*k-16)= al(36*k-16)/(scale(ip4)*scale(iq2))
362  al(36*k-15)= al(36*k-15)/(scale(ip4)*scale(iq3))
363  al(36*k-14)= al(36*k-14)/(scale(ip4)*scale(iq4))
364  al(36*k-13)= al(36*k-13)/(scale(ip4)*scale(iq5))
365  al(36*k-12)= al(36*k-12)/(scale(ip4)*scale(iq6))
366 
367  al(36*k-11)= al(36*k-11)/(scale(ip5)*scale(iq1))
368  al(36*k-10)= al(36*k-10)/(scale(ip5)*scale(iq2))
369  al(36*k-9 )= al(36*k-9 )/(scale(ip5)*scale(iq3))
370  al(36*k-8 )= al(36*k-8 )/(scale(ip5)*scale(iq4))
371  al(36*k-7 )= al(36*k-7 )/(scale(ip5)*scale(iq5))
372  al(36*k-6 )= al(36*k-6 )/(scale(ip5)*scale(iq6))
373 
374  al(36*k-5 )= al(36*k-5 )/(scale(ip6)*scale(iq1))
375  al(36*k-4 )= al(36*k-4 )/(scale(ip6)*scale(iq2))
376  al(36*k-3 )= al(36*k-3 )/(scale(ip6)*scale(iq3))
377  al(36*k-2 )= al(36*k-2 )/(scale(ip6)*scale(iq4))
378  al(36*k-1 )= al(36*k-1 )/(scale(ip6)*scale(iq5))
379  al(36*k )= al(36*k )/(scale(ip6)*scale(iq6))
380  enddo
381 
382  isu= inu(i-1) + 1
383  ieu= inu(i )
384  !*voption indep (IAU,AU,SCALE)
385  do k= isu, ieu
386  inod= iau(k)
387  iq1= 6*inod - 5
388  iq2= 6*inod - 4
389  iq3= 6*inod - 3
390  iq4= 6*inod - 2
391  iq5= 6*inod - 1
392  iq6= 6*inod
393 
394  au(36*k-35)= au(36*k-35)/(scale(ip1)*scale(iq1))
395  au(36*k-34)= au(36*k-34)/(scale(ip1)*scale(iq2))
396  au(36*k-33)= au(36*k-33)/(scale(ip1)*scale(iq3))
397  au(36*k-32)= au(36*k-32)/(scale(ip1)*scale(iq4))
398  au(36*k-31)= au(36*k-31)/(scale(ip1)*scale(iq5))
399  au(36*k-30)= au(36*k-30)/(scale(ip1)*scale(iq6))
400 
401  au(36*k-29)= au(36*k-29)/(scale(ip2)*scale(iq1))
402  au(36*k-28)= au(36*k-28)/(scale(ip2)*scale(iq2))
403  au(36*k-27)= au(36*k-27)/(scale(ip2)*scale(iq3))
404  au(36*k-26)= au(36*k-26)/(scale(ip2)*scale(iq4))
405  au(36*k-25)= au(36*k-25)/(scale(ip2)*scale(iq5))
406  au(36*k-24)= au(36*k-24)/(scale(ip2)*scale(iq6))
407 
408  au(36*k-23)= au(36*k-23)/(scale(ip3)*scale(iq1))
409  au(36*k-22)= au(36*k-22)/(scale(ip3)*scale(iq2))
410  au(36*k-21)= au(36*k-21)/(scale(ip3)*scale(iq3))
411  au(36*k-20)= au(36*k-20)/(scale(ip3)*scale(iq4))
412  au(36*k-19)= au(36*k-19)/(scale(ip3)*scale(iq5))
413  au(36*k-18)= au(36*k-18)/(scale(ip3)*scale(iq6))
414 
415  au(36*k-17)= au(36*k-17)/(scale(ip4)*scale(iq1))
416  au(36*k-16)= au(36*k-16)/(scale(ip4)*scale(iq2))
417  au(36*k-15)= au(36*k-15)/(scale(ip4)*scale(iq3))
418  au(36*k-14)= au(36*k-14)/(scale(ip4)*scale(iq4))
419  au(36*k-13)= au(36*k-13)/(scale(ip4)*scale(iq5))
420  au(36*k-12)= au(36*k-12)/(scale(ip4)*scale(iq6))
421 
422  au(36*k-11)= au(36*k-11)/(scale(ip5)*scale(iq1))
423  au(36*k-10)= au(36*k-10)/(scale(ip5)*scale(iq2))
424  au(36*k-9 )= au(36*k-9 )/(scale(ip5)*scale(iq3))
425  au(36*k-8 )= au(36*k-8 )/(scale(ip5)*scale(iq4))
426  au(36*k-7 )= au(36*k-7 )/(scale(ip5)*scale(iq5))
427  au(36*k-6 )= au(36*k-6 )/(scale(ip5)*scale(iq6))
428 
429  au(36*k-5 )= au(36*k-5 )/(scale(ip6)*scale(iq1))
430  au(36*k-4 )= au(36*k-4 )/(scale(ip6)*scale(iq2))
431  au(36*k-3 )= au(36*k-3 )/(scale(ip6)*scale(iq3))
432  au(36*k-2 )= au(36*k-2 )/(scale(ip6)*scale(iq4))
433  au(36*k-1 )= au(36*k-1 )/(scale(ip6)*scale(iq5))
434  au(36*k )= au(36*k )/(scale(ip6)*scale(iq6))
435  enddo
436  enddo
437 
438  deallocate(scale)
439  end subroutine hecmw_solver_scaling_bk_66
440 
441 end module hecmw_solver_scaling_66
hecmw_util::hecmw_wtime
real(kind=kreal) function hecmw_wtime()
Definition: hecmw_util_f.F90:549
hecmw_solver_scaling_66
Definition: hecmw_solver_scaling_66.f90:6
hecmw_solver_scaling_66::hecmw_solver_scaling_bk_66
subroutine, public hecmw_solver_scaling_bk_66(hecMAT)
Definition: hecmw_solver_scaling_66.f90:235
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
hecmw_util::hecmwst_local_mesh
Definition: hecmw_util_f.F90:234
m_hecmw_comm_f
Definition: hecmw_comm_f.F90:6
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
hecmw_solver_scaling_66::hecmw_solver_scaling_fw_66
subroutine, public hecmw_solver_scaling_fw_66(hecMESH, hecMAT, COMMtime)
Definition: hecmw_solver_scaling_66.f90:21
hecmw_matrix_misc
Definition: hecmw_matrix_misc.f90:6
hecmw_matrix_misc::hecmw_mat_get_scaling
integer(kind=kint) function, public hecmw_mat_get_scaling(hecMAT)
Definition: hecmw_matrix_misc.f90:394
m_hecmw_comm_f::hecmw_update_r
subroutine hecmw_update_r(hecMESH, val, n, m)
Definition: hecmw_comm_f.F90:683
hecmw_util::hecmwst_matrix
Definition: hecmw_util_f.F90:444