FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hecmw_mat_ass.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 
14  public :: hecmw_mat_ass_elem
15  public :: hecmw_mat_add_node
16  public :: hecmw_array_search_i
17  public :: hecmw_mat_ass_equation
19  public :: hecmw_mat_add_dof
20  public :: hecmw_mat_ass_bc
22  public :: hecmw_mat_ass_contactlag
23  public :: stf_get_block
24 
25 contains
26 
27  !C
28  !C***
29  !C*** MAT_ASS_ELEM
30  !C***
31  !C
32  subroutine hecmw_mat_ass_elem(hecMAT, nn, nodLOCAL, stiffness)
33  type (hecmwst_matrix) :: hecmat
34  integer(kind=kint) :: nn
35  integer(kind=kint) :: nodlocal(:)
36  real(kind=kreal) :: stiffness(:, :)
37  !** Local variables
38  integer(kind=kint) :: ndof, inod_e, jnod_e, inod, jnod
39  real(kind=kreal) :: a(6,6)
40 
41  ndof = hecmat%NDOF
42 
43  do inod_e = 1, nn
44  inod = nodlocal(inod_e)
45  do jnod_e = 1, nn
46  jnod = nodlocal(jnod_e)
47  !***** Add components
48  call stf_get_block(stiffness, ndof, inod_e, jnod_e, a)
49  call hecmw_mat_add_node(hecmat, inod, jnod, a)
50  enddo
51  enddo
52 
53  end subroutine hecmw_mat_ass_elem
54 
55  subroutine stf_get_block(stiffness, ndof, inod, jnod, a)
56  real(kind=kreal) :: stiffness(:, :), a(:, :)
57  integer(kind=kint) :: ndof, inod, jnod
58  !** Local variables
59  integer(kind=kint) :: row_offset, col_offset, i, j
60 
61  row_offset = ndof*(inod-1)
62  do i = 1, ndof
63 
64  col_offset = ndof*(jnod-1)
65  do j = 1, ndof
66 
67  a(i, j) = stiffness(i + row_offset, j + col_offset)
68  enddo
69  enddo
70  end subroutine stf_get_block
71 
72 
73  subroutine hecmw_mat_add_node(hecMAT, inod, jnod, a)
74  type (hecmwst_matrix) :: hecmat
75  integer(kind=kint) :: inod, jnod
76  real(kind=kreal) :: a(:, :)
77  !** Local variables
78  integer(kind=kint) :: ndof, is, ie, k, idx_base, idx, idof, jdof
79 
80  ndof = hecmat%NDOF
81 
82  if (inod < jnod) then
83  is = hecmat%indexU(inod-1)+1
84  ie = hecmat%indexU(inod)
85  k = hecmw_array_search_i(hecmat%itemU, is, ie, jnod)
86 
87  if (k < is .or. ie < k) then
88  write(*,*) '###ERROR### : cannot find connectivity (1)'
89  write(*,*) ' myrank = ', hecmw_comm_get_rank(), ', inod = ', inod, ', jnod = ', jnod
91  endif
92 
93  idx_base = ndof**2 * (k-1)
94  do idof = 1, ndof
95  do jdof = 1, ndof
96  idx = idx_base + jdof
97  !$omp atomic
98  hecmat%AU(idx) = hecmat%AU(idx) + a(idof, jdof)
99  enddo
100  idx_base = idx_base + ndof
101  enddo
102 
103  else if (inod > jnod) then
104  is = hecmat%indexL(inod-1)+1
105  ie = hecmat%indexL(inod)
106  k = hecmw_array_search_i(hecmat%itemL, is, ie, jnod)
107 
108  if (k < is .or. ie < k) then
109  write(*,*) '###ERROR### : cannot find connectivity (2)'
110  write(*,*) ' myrank = ', hecmw_comm_get_rank(), ', inod = ', inod, ', jnod = ', jnod
112  endif
113 
114  idx_base = ndof**2 * (k-1)
115  do idof = 1, ndof
116  do jdof = 1, ndof
117  idx = idx_base + jdof
118  !$omp atomic
119  hecmat%AL(idx) = hecmat%AL(idx) + a(idof, jdof)
120  enddo
121  idx_base = idx_base + ndof
122  enddo
123 
124  else
125  idx_base = ndof**2 * (inod - 1)
126  do idof = 1, ndof
127  do jdof = 1, ndof
128  idx = idx_base + jdof
129  !$omp atomic
130  hecmat%D(idx) = hecmat%D(idx) + a(idof, jdof)
131  enddo
132  idx_base = idx_base + ndof
133  enddo
134  endif
135  end subroutine hecmw_mat_add_node
136 
137 
138  function hecmw_array_search_i(array, is, iE, ival)
139  integer(kind=kint) :: hecmw_array_search_i
140  integer(kind=kint) :: array(:)
141  integer(kind=kint) :: is, ie, ival
142  !** Local variables
143  integer(kind=kint) :: left, right, center, cval
144 
145  left = is
146  right = ie
147  do
148  if (left > right) then
149  center = -1
150  exit
151  endif
152 
153  center = (left + right) / 2
154  cval = array(center)
155 
156  if (ival < cval) then
157  right = center - 1
158  else if (cval < ival) then
159  left = center + 1
160  else
161  exit
162  endif
163  enddo
164 
165  hecmw_array_search_i = center
166 
167  end function hecmw_array_search_i
168 
169 
170  !C
171  !C***
172  !C*** MAT_ASS_EQUATION
173  !C***
174  !C
175  subroutine hecmw_mat_ass_equation ( hecMESH, hecMAT )
176  type (hecmwst_matrix), target :: hecmat
177  type (hecmwst_local_mesh) :: hecmesh
178  !** Local variables
179  real(kind=kreal), pointer :: penalty
180  real(kind=kreal) :: alpha, a1_2inv, ai, aj, factor
181  integer(kind=kint) :: impc, is, ie, i, j, inod, idof, jnod, jdof
182  logical :: is_internal_i
183 
184  if( hecmw_mat_get_penalized(hecmat) == 1 ) return
185 
186  ! write(*,*) "INFO: imposing MPC by penalty"
187 
188  penalty => hecmat%Rarray(11)
189 
190  if (penalty < 0.0) stop "ERROR: negative penalty"
191  if (penalty < 1.0) write(*,*) "WARNING: penalty ", penalty, " smaller than 1"
192 
193  alpha= hecmw_mat_diag_max(hecmat, hecmesh) * penalty
194  call hecmw_mat_set_penalty_alpha(hecmat, alpha)
195 
196  outer: do impc = 1, hecmesh%mpc%n_mpc
197  is = hecmesh%mpc%mpc_index(impc-1) + 1
198  ie = hecmesh%mpc%mpc_index(impc)
199 
200  do i = is, ie
201  if (hecmesh%mpc%mpc_dof(i) > hecmat%NDOF) cycle outer
202  enddo
203 
204  a1_2inv = 1.0 / hecmesh%mpc%mpc_val(is)**2
205 
206 
207  do i = is, ie
208  inod = hecmesh%mpc%mpc_item(i)
209 
210  is_internal_i = (hecmesh%node_ID(2*inod) == hecmw_comm_get_rank())
211  if (.not. is_internal_i) cycle
212 
213  idof = hecmesh%mpc%mpc_dof(i)
214  ai = hecmesh%mpc%mpc_val(i)
215  factor = ai * a1_2inv
216 
217  do j = is, ie
218  jnod = hecmesh%mpc%mpc_item(j)
219  jdof = hecmesh%mpc%mpc_dof(j)
220  aj = hecmesh%mpc%mpc_val(j)
221 
222  call hecmw_mat_add_dof(hecmat, inod, idof, jnod, jdof, aj*factor*alpha)
223  enddo
224 
225  enddo
226  enddo outer
227 
228  call hecmw_mat_set_penalized(hecmat, 1)
229 
230  end subroutine hecmw_mat_ass_equation
231 
232 
233  subroutine hecmw_mat_ass_equation_rhs ( hecMESH, hecMAT )
234  type (hecmwst_matrix), target :: hecmat
235  type (hecmwst_local_mesh) :: hecmesh
236  !** Local variables
237  real(kind=kreal) :: alpha, a1_2inv, ai, factor, ci
238  integer(kind=kint) :: ndof, impc, is, ie, i, inod, idof
239 
240  if( hecmw_mat_get_penalized_b(hecmat) == 1) return
241 
242  alpha = hecmw_mat_get_penalty_alpha(hecmat)
243  if (alpha <= 0.0) stop "ERROR: penalty applied on vector before matrix"
244 
245  ndof = hecmat%NDOF
246 
247  outer: do impc = 1, hecmesh%mpc%n_mpc
248  is = hecmesh%mpc%mpc_index(impc-1) + 1
249  ie = hecmesh%mpc%mpc_index(impc)
250 
251  do i = is, ie
252  if (hecmesh%mpc%mpc_dof(i) > ndof) cycle outer
253  enddo
254 
255  a1_2inv = 1.0 / hecmesh%mpc%mpc_val(is)**2
256 
257  do i = is, ie
258  inod = hecmesh%mpc%mpc_item(i)
259 
260  idof = hecmesh%mpc%mpc_dof(i)
261  ai = hecmesh%mpc%mpc_val(i)
262  factor = ai * a1_2inv
263 
264  ci = hecmesh%mpc%mpc_const(impc)
265  !$omp atomic
266  hecmat%B(ndof*(inod-1)+idof) = hecmat%B(ndof*(inod-1)+idof) + ci*factor*alpha
267  enddo
268  enddo outer
269 
270  call hecmw_mat_set_penalized_b(hecmat, 1)
271 
272  end subroutine hecmw_mat_ass_equation_rhs
273 
274 
275  subroutine hecmw_mat_add_dof(hecMAT, inod, idof, jnod, jdof, val)
276  type (hecmwst_matrix) :: hecmat
277  integer(kind=kint) :: inod, idof, jnod, jdof
278  real(kind=kreal) :: val
279  !** Local variables
280  integer(kind=kint) :: ndof, is, ie, k, idx
281 
282  ndof = hecmat%NDOF
283  if (inod < jnod) then
284  is = hecmat%indexU(inod-1)+1
285  ie = hecmat%indexU(inod)
286  k = hecmw_array_search_i(hecmat%itemU, is, ie, jnod)
287 
288  if (k < is .or. ie < k) then
289  write(*,*) '###ERROR### : cannot find connectivity (3)'
290  write(*,*) ' myrank = ', hecmw_comm_get_rank(), ', inod = ', inod, ', jnod = ', jnod
292  return
293  endif
294 
295  idx = ndof**2 * (k-1) + ndof * (idof-1) + jdof
296  !$omp atomic
297  hecmat%AU(idx) = hecmat%AU(idx) + val
298 
299  else if (inod > jnod) then
300  is = hecmat%indexL(inod-1)+1
301  ie = hecmat%indexL(inod)
302  k = hecmw_array_search_i(hecmat%itemL, is, ie, jnod)
303 
304  if (k < is .or. ie < k) then
305  write(*,*) '###ERROR### : cannot find connectivity (4)'
306  write(*,*) ' myrank = ', hecmw_comm_get_rank(), ', inod = ', inod, ', jnod = ', jnod
308  return
309  endif
310 
311  idx = ndof**2 * (k-1) + ndof * (idof-1) + jdof
312  !$omp atomic
313  hecmat%AL(idx) = hecmat%AL(idx) + val
314 
315  else
316  idx = ndof**2 * (inod - 1) + ndof * (idof - 1) + jdof
317  !$omp atomic
318  hecmat%D(idx) = hecmat%D(idx) + val
319  endif
320 
321  end subroutine hecmw_mat_add_dof
322 
323  !C
324  !C***
325  !C*** MAT_ASS_BC
326  !C***
327  !C
328  subroutine hecmw_mat_ass_bc(hecMAT, inode, idof, RHS, conMAT)
329  type (hecmwst_matrix) :: hecmat
330  integer(kind=kint) :: inode, idof
331  real(kind=kreal) :: rhs, val
332  type (hecmwst_matrix),optional :: conmat
333  integer(kind=kint) :: ndof, in, i, ii, iii, ndof2, k, is, ie, iis, iie, ik, idx
334 
335  ndof = hecmat%NDOF
336  if( ndof < idof ) return
337 
338  !C-- DIAGONAL block
339 
340  hecmat%B(ndof*inode-(ndof-idof)) = rhs
341  if(present(conmat)) conmat%B(ndof*inode-(ndof-idof)) = 0.0d0
342  ndof2 = ndof*ndof
343  ii = ndof2 - idof
344 
345  do i = ndof-1,0,-1
346  if( i .NE. ndof-idof ) then
347  idx = ndof*inode-i
348  val = hecmat%D(ndof2*inode-ii)*rhs
349  !$omp atomic
350  hecmat%B(idx) = hecmat%B(idx) - val
351  if(present(conmat)) then
352  val = conmat%D(ndof2*inode-ii)*rhs
353  !$omp atomic
354  conmat%B(idx) = conmat%B(idx) - val
355  endif
356  endif
357  ii = ii - ndof
358  end do
359 
360  !*Set diagonal row to zero
361  ii = ndof2-1 - (idof-1)*ndof
362 
363  do i = 0, ndof - 1
364  hecmat%D(ndof2*inode-ii+i)=0.d0
365  if(present(conmat)) conmat%D(ndof2*inode-ii+i)=0.d0
366  end do
367 
368  !*Set diagonal column to zero
369  ii = ndof2 - idof
370  do i = 1, ndof
371  if( i.NE.idof ) then
372  hecmat%D(ndof2*inode-ii) = 0.d0
373  if(present(conmat)) conmat%D(ndof2*inode-ii) = 0.d0
374  else
375  hecmat%D(ndof2*inode-ii) = 1.d0
376  if(present(conmat)) conmat%D(ndof2*inode-ii) = 0.d0
377  endif
378  ii = ii - ndof
379  end do
380 
381  !C-- OFF-DIAGONAL blocks
382 
383  ii = ndof2-1 - (idof-1)*ndof
384  is = hecmat%indexL(inode-1) + 1
385  ie = hecmat%indexL(inode )
386 
387  do k= is, ie
388 
389  !*row (left)
390  do i = 0, ndof - 1
391  hecmat%AL(ndof2*k-ii+i) = 0.d0
392  if(present(conmat)) conmat%AL(ndof2*k-ii+i) = 0.d0
393  end do
394 
395  !*column (upper)
396  in = hecmat%itemL(k)
397  iis = hecmat%indexU(in-1) + 1
398  iie = hecmat%indexU(in )
399  do ik = iis, iie
400  if (hecmat%itemU(ik) .eq. inode) then
401  iii = ndof2 - idof
402  do i = ndof-1,0,-1
403  idx = ndof*in-i
404  val = hecmat%AU(ndof2*ik-iii)*rhs
405  !$omp atomic
406  hecmat%B(idx) = hecmat%B(idx) - val
407  hecmat%AU(ndof2*ik-iii)= 0.d0
408  if(present(conmat)) then
409  val = conmat%AU(ndof2*ik-iii)*rhs
410  !$omp atomic
411  conmat%B(idx) = conmat%B(idx) - val
412  conmat%AU(ndof2*ik-iii)= 0.d0
413  endif
414  iii = iii - ndof
415  end do
416  exit
417  endif
418  enddo
419 
420  enddo
421 
422  ii = ndof2-1 - (idof-1)*ndof
423  is = hecmat%indexU(inode-1) + 1
424  ie = hecmat%indexU(inode )
425 
426  do k= is, ie
427 
428  !*row (right)
429  do i = 0,ndof-1
430  hecmat%AU(ndof2*k-ii+i) = 0.d0
431  if(present(conmat)) conmat%AU(ndof2*k-ii+i) = 0.d0
432  end do
433 
434  !*column (lower)
435  in = hecmat%itemU(k)
436  iis = hecmat%indexL(in-1) + 1
437  iie = hecmat%indexL(in )
438  do ik= iis, iie
439  if (hecmat%itemL(ik) .eq. inode) then
440  iii = ndof2 - idof
441 
442  do i = ndof-1, 0, -1
443  idx = ndof*in-i
444  val = hecmat%AL(ndof2*ik-iii)*rhs
445  !$omp atomic
446  hecmat%B(idx) = hecmat%B(idx) - val
447  hecmat%AL(ndof2*ik-iii) = 0.d0
448  if(present(conmat)) then
449  val = conmat%AL(ndof2*ik-iii)*rhs
450  !$omp atomic
451  conmat%B(idx) = conmat%B(idx) - val
452  conmat%AL(ndof2*ik-iii) = 0.d0
453  endif
454  iii = iii - ndof
455  end do
456  exit
457  endif
458  enddo
459 
460  enddo
461  !*End off - diagonal blocks
462 
463  end subroutine hecmw_mat_ass_bc
464 
465 
468  subroutine hecmw_mat_ass_bc_contactlag(hecMAT,hecLagMAT,inode,idof,RHS)
469 
470  type(hecmwst_matrix) :: hecmat
471  type(hecmwst_matrix_lagrange) :: heclagmat
472  integer(kind=kint) :: inode, idof
473  integer(kind=kint) :: isu, ieu, isl, iel, i, l, k
474  real(kind=kreal) :: rhs
475 
476  isu = heclagmat%indexU_lagrange(inode-1)+1
477  ieu = heclagmat%indexU_lagrange(inode)
478  do i = isu, ieu
479  heclagmat%AU_lagrange((i-1)*3+idof) = 0.0d0
480  l = heclagmat%itemU_lagrange(i)
481  isl = heclagmat%indexL_lagrange(l-1)+1
482  iel = heclagmat%indexL_lagrange(l)
483  k = hecmw_array_search_i(heclagmat%itemL_lagrange,isl,iel,inode)
484  if(k < isl .or. k > iel) cycle
485  hecmat%B(hecmat%NP*hecmat%NDOF+l) = hecmat%B(hecmat%NP*hecmat%NDOF+l) - heclagmat%AL_lagrange((k-1)*3+idof)*rhs
486  heclagmat%AL_lagrange((k-1)*3+idof) = 0.0d0
487  enddo
488 
489  end subroutine hecmw_mat_ass_bc_contactlag
490 
492  subroutine hecmw_mat_ass_contactlag(nnode,ndLocal,id_lagrange,fcoeff,stiffness,hecMAT,hecLagMAT)
493 
494  type(hecmwst_matrix) :: hecmat
495  type(hecmwst_matrix_lagrange) :: heclagmat
496  integer(kind=kint) :: nnode, ndlocal(nnode + 1), id_lagrange
499  integer(kind=kint) :: i, j, inod, jnod, l
500  integer(kind=kint) :: isl, iel, idxl_base, kl, idxl, isu, ieu, idxu_base, ku, idxu
501  real(kind=kreal) :: fcoeff
502  real(kind=kreal) :: stiffness(:,:)
503  real(kind=kreal) :: a(3, 3)
504 
505  i = nnode + 1 + 1
506  inod = id_lagrange
507  isl = heclagmat%indexL_lagrange(inod-1)+1
508  iel = heclagmat%indexL_lagrange(inod)
509 
510  do j = 1, nnode + 1
511  jnod = ndlocal(j)
512  isu = heclagmat%indexU_lagrange(jnod-1)+1
513  ieu = heclagmat%indexU_lagrange(jnod)
514 
515  kl = hecmw_array_search_i(heclagmat%itemL_lagrange,isl,iel,jnod)
516  if( kl<isl .or. kl>iel ) then
517  write(*,*) '###ERROR### : cannot find connectivity (Lagrange1)'
518  stop
519  endif
520  ku = hecmw_array_search_i(heclagmat%itemU_lagrange,isu,ieu,inod)
521  if( ku<isu .or. ku>ieu ) then
522  write(*,*) '###ERROR### : cannot find connectivity (Lagrange2)'
523  stop
524  endif
525 
526  idxl_base = (kl-1)*3
527  idxu_base = (ku-1)*3
528 
529  do l = 1, 3
530  idxl = idxl_base + l
531  heclagmat%AL_lagrange(idxl) = heclagmat%AL_lagrange(idxl) + stiffness((i-1)*3+1,(j-1)*3+l)
532  idxu = idxu_base + l
533  heclagmat%AU_lagrange(idxu) = heclagmat%AU_lagrange(idxu) + stiffness((j-1)*3+l,(i-1)*3+1)
534  enddo
535  enddo
536 
537 
538  if(fcoeff /= 0.0d0)then
539 
540  do i = 1, nnode + 1
541  inod = ndlocal(i)
542  do j = 1, nnode + 1
543  jnod = ndlocal(j)
544  call stf_get_block(stiffness(1:(nnode+1)*3,1:(nnode+1)*3), 3, i, j, a)
545  call hecmw_mat_add_node(hecmat, inod, jnod, a)
546  enddo
547  enddo
548 
549  endif
550 
551  end subroutine hecmw_mat_ass_contactlag
552 
553 end module hecmw_matrix_ass
hecmw_matrix_misc::hecmw_mat_get_penalty_alpha
real(kind=kreal) function, public hecmw_mat_get_penalty_alpha(hecMAT)
Definition: hecmw_matrix_misc.f90:770
hecmw_util::hecmw_abort
subroutine hecmw_abort(comm)
Definition: hecmw_util_f.F90:534
hecmw_matrix_ass::hecmw_mat_ass_elem
subroutine, public hecmw_mat_ass_elem(hecMAT, nn, nodLOCAL, stiffness)
Definition: hecmw_mat_ass.f90:33
hecmw_matrix_ass::hecmw_mat_ass_bc_contactlag
subroutine, public hecmw_mat_ass_bc_contactlag(hecMAT, hecLagMAT, inode, idof, RHS)
Modify Lagrange multiplier-related part of stiffness matrix and right-hand side vector for dealing wi...
Definition: hecmw_mat_ass.f90:469
hecmw_matrix_misc::hecmw_mat_diag_max
real(kind=kreal) function, public hecmw_mat_diag_max(hecMAT, hecMESH)
Definition: hecmw_matrix_misc.f90:777
hecmw_matrix_misc::hecmw_mat_set_penalty_alpha
subroutine, public hecmw_mat_set_penalty_alpha(hecMAT, alpha)
Definition: hecmw_matrix_misc.f90:763
hecmw_matrix_ass::hecmw_array_search_i
integer(kind=kint) function, public hecmw_array_search_i(array, is, iE, ival)
Definition: hecmw_mat_ass.f90:139
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
hecmw_util::hecmwst_local_mesh
Definition: hecmw_util_f.F90:234
hecmw_matrix_misc::hecmw_mat_set_penalized
subroutine, public hecmw_mat_set_penalized(hecMAT, penalized)
Definition: hecmw_matrix_misc.f90:401
hecmw_matrix_misc::hecmw_mat_get_penalized
integer(kind=kint) function, public hecmw_mat_get_penalized(hecMAT)
Definition: hecmw_matrix_misc.f90:408
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_matrix_ass::hecmw_mat_ass_bc
subroutine, public hecmw_mat_ass_bc(hecMAT, inode, idof, RHS, conMAT)
Definition: hecmw_mat_ass.f90:329
hecmw_matrix_misc::hecmw_mat_set_penalized_b
subroutine, public hecmw_mat_set_penalized_b(hecMAT, penalized_b)
Definition: hecmw_matrix_misc.f90:415
hecmw_matrix_ass
Definition: hecmw_mat_ass.f90:6
hecmw_matrix_misc::hecmw_mat_get_penalized_b
integer(kind=kint) function, public hecmw_mat_get_penalized_b(hecMAT)
Definition: hecmw_matrix_misc.f90:422
hecmw_matrix_misc
Definition: hecmw_matrix_misc.f90:6
hecmw_matrix_ass::hecmw_mat_ass_contactlag
subroutine, public hecmw_mat_ass_contactlag(nnode, ndLocal, id_lagrange, fcoeff, stiffness, hecMAT, hecLagMAT)
This subroutine assembles contact stiffness matrix of a contact pair into global stiffness matrix.
Definition: hecmw_mat_ass.f90:493
hecmw_matrix_ass::stf_get_block
subroutine, public stf_get_block(stiffness, ndof, inod, jnod, a)
Definition: hecmw_mat_ass.f90:56
hecmw_matrix_ass::hecmw_mat_ass_equation
subroutine, public hecmw_mat_ass_equation(hecMESH, hecMAT)
Definition: hecmw_mat_ass.f90:176
hecmw_matrix_ass::hecmw_mat_add_node
subroutine, public hecmw_mat_add_node(hecMAT, inod, jnod, a)
Definition: hecmw_mat_ass.f90:74
hecmw_util::hecmwst_matrix_lagrange
Structure for Lagrange multiplier-related part of stiffness matrix (Lagrange multiplier-related matri...
Definition: hecmw_util_f.F90:427
hecmw_util::hecmw_comm_get_rank
integer(kind=kint) function hecmw_comm_get_rank()
Definition: hecmw_util_f.F90:582
hecmw_matrix_ass::hecmw_mat_ass_equation_rhs
subroutine, public hecmw_mat_ass_equation_rhs(hecMESH, hecMAT)
Definition: hecmw_mat_ass.f90:234
hecmw_util::hecmw_comm_get_comm
integer(kind=kint) function hecmw_comm_get_comm()
Definition: hecmw_util_f.F90:571
hecmw_matrix_ass::hecmw_mat_add_dof
subroutine, public hecmw_mat_add_dof(hecMAT, inod, idof, jnod, jdof, val)
Definition: hecmw_mat_ass.f90:276
hecmw_util::hecmwst_matrix
Definition: hecmw_util_f.F90:444