FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hecmw_solver_contact.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 !-------------------------------------------------------------------------------
8 
10 
11  use hecmw_util
18  use m_hecmw_comm_f
19 
20  implicit none
21 
22  private
23  public :: solve_lineq_contact_init
24  public :: solve_lineq_contact
25 
26 contains
27 
29  subroutine solve_lineq_contact_init(hecMESH,hecMAT,hecLagMAT,is_sym)
30  type (hecmwst_local_mesh) :: hecmesh
31  type (hecmwst_matrix) :: hecmat
32  type (hecmwst_matrix_lagrange) :: heclagmat
33  logical :: is_sym
34  integer(kind=kint) :: contact_elim, solver_type, iterlog, timelog, myrank
35  logical :: writelog
36 
37  contact_elim = hecmw_mat_get_contact_elim(hecmat)
38  solver_type = hecmw_mat_get_solver_type(hecmat)
39  iterlog = hecmw_mat_get_iterlog(hecmat)
40  timelog = hecmw_mat_get_timelog(hecmat)
42  writelog = (myrank==0 .and. (iterlog>0 .or. timelog>0))
43 
44  if( contact_elim==0 )then ! auto
45  if( solver_type==1 )then ! iterative
46  contact_elim=1
47  else ! direct
48  contact_elim=-1
49  endif
50  call hecmw_mat_set_contact_elim(hecmat,contact_elim)
51  endif
52 
53  if( contact_elim==1 )then
54  if( writelog ) write(*,*) 'solve contact with elimination'
55  call solve_lineq_contact_elim_init(hecmesh,hecmat,heclagmat,is_sym)
56  else
57  if( writelog ) write(*,*) 'solve contact without elimination'
58  if( solver_type==1 )then
59  write(*,*) 'ERROR: iterative solver without elimination not available in contact analysis'
61  elseif( solver_type==2 )then
62  call solve_lineq_serial_lag_hecmw_init(hecmat,heclagmat,is_sym)
63  else if( solver_type==3 )then
64  call solve_lineq_mkl_contact_init(hecmesh,is_sym)
65  elseif( solver_type==5 ) then
66  call solve_lineq_mumps_contact_init(hecmesh,hecmat,heclagmat,is_sym)
67  endif
68  endif
69  end subroutine solve_lineq_contact_init
70 
71 
73  subroutine solve_lineq_contact(hecMESH,hecMAT,hecLagMAT,conMAT,istat,rf,is_contact_active)
74 
75  type (hecmwst_local_mesh) :: hecmesh
76  type (hecmwst_matrix) :: hecmat
77  type (hecmwst_matrix_lagrange) :: heclagmat
78  type (hecmwst_matrix) :: conmat
79  integer(kind=kint), intent(out) :: istat
80  real(kind=kreal), optional :: rf
81  logical :: is_contact_active
82 
83  real(kind=kreal) :: factor
84  real(kind=kreal) :: t1, t2
85  integer(kind=kint) :: ndof
86  integer(kind=kint) :: contact_elim, solver_type
87 
88  contact_elim = hecmw_mat_get_contact_elim(hecmat)
89  solver_type = hecmw_mat_get_solver_type(hecmat)
90 
91  factor = 1.0d0
92  if( present(rf) )factor = rf
93 
94  t1 = hecmw_wtime()
95 
96  istat = 0
97  if( contact_elim==1 )then
98  call solve_lineq_contact_elim(hecmesh,hecmat,heclagmat,istat,conmat,is_contact_active)
99  else
100  if( solver_type==1 )then
101  write(*,*) 'ERROR: iterative solver without elimination not available in contact analysis'
103  elseif( solver_type==2 )then
104  if( hecmw_comm_get_size() > 1) then
105  write(*,*) 'ERROR: !SOLVER,METHOD=DIRECT not available in parallel contact analysis;',&
106  ' please use MUMPS or DIRECTmkl instead'
108  else
109  call add_conmat_to_hecmat(hecmat,conmat,heclagmat)
110  call solve_lineq_serial_lag_hecmw(hecmesh,hecmat,heclagmat)
111  endif
112  elseif( solver_type==3 )then
113  if( hecmw_comm_get_size() > 1) then
114  call solve_lineq_mkl_contact(hecmesh,hecmat,heclagmat,istat,conmat)
115  else
116  call add_conmat_to_hecmat(hecmat,conmat,heclagmat)
117  call solve_lineq_mkl_contact(hecmesh,hecmat,heclagmat,istat)
118  endif
119  elseif( solver_type==5 ) then
120  call solve_lineq_mumps_contact(hecmesh,hecmat,heclagmat,istat,conmat)
121  endif
122  endif
123 
124  ndof = hecmat%NDOF
125  call hecmw_update_r(hecmesh,hecmat%X,hecmesh%n_node, ndof)
126 
127  t2 = hecmw_wtime()
128  if (hecmw_mat_get_timelog(hecmat) .ge. 1) then
129  if ( hecmw_comm_get_rank() ==0) write(*,*) ' solve time :', t2 - t1
130  endif
131 
132  hecmat%X=factor*hecmat%X
133 
134  end subroutine solve_lineq_contact
135 
136 
137  subroutine add_conmat_to_hecmat(hecMAT,conMAT,hecLagMat)
138  type (hecmwst_matrix), intent(inout) :: hecmat
139  type (hecmwst_matrix), intent(in) :: conmat
140  type (hecmwst_matrix_lagrange), intent(in) :: heclagmat
141 
142 
143  integer(kind=kint) :: ndof,ndof2,i
144 
145  ndof = hecmat%NDOF
146  ndof2 = ndof*ndof
147 
148  do i=1,hecmat%NP*ndof + heclagmat%num_lagrange
149  hecmat%B(i) = hecmat%B(i) + conmat%B(i)
150  enddo
151 
152  do i=1,hecmat%NP*ndof2
153  hecmat%D(i) = hecmat%D(i) + conmat%D(i)
154  enddo
155 
156  do i=1,hecmat%NPL*ndof2
157  hecmat%AL(i) = hecmat%AL(i) + conmat%AL(i)
158  enddo
159 
160  do i=1,hecmat%NPU*ndof2
161  hecmat%AU(i) = hecmat%AU(i) + conmat%AU(i)
162  enddo
163  end subroutine add_conmat_to_hecmat
164 
165 
166 end module m_solve_lineq_contact
hecmw_matrix_misc::hecmw_mat_get_solver_type
integer(kind=kint) function, public hecmw_mat_get_solver_type(hecMAT)
Definition: hecmw_matrix_misc.f90:601
m_solve_lineq_mumps_contact::solve_lineq_mumps_contact
subroutine, public solve_lineq_mumps_contact(hecMESH, hecMAT, hecLagMAT, istat, conMAT)
Definition: solve_LINEQ_MUMPS_contact.f90:68
m_solve_lineq_contact_elim
This module provides interface of iteratie linear equation solver for contact problems using Lagrange...
Definition: solve_LINEQ_contact_elim.f90:7
m_hecmw_mat_resid_contact
Definition: hecmw_mat_resid_contact.f90:5
hecmw_util::hecmw_wtime
real(kind=kreal) function hecmw_wtime()
Definition: hecmw_util_f.F90:549
m_solve_lineq_mumps_contact
This module provides linear equation solver interface of MUMPS for contact problems using Lagrange mu...
Definition: solve_LINEQ_MUMPS_contact.f90:7
hecmw_util::hecmw_abort
subroutine hecmw_abort(comm)
Definition: hecmw_util_f.F90:534
m_solve_lineq_mkl_contact
This module provides functions to solve sparse system of \linear equitions using intel MKL direct spa...
Definition: solve_LINEQ_mkl.F90:8
hecmw_matrix_misc::hecmw_mat_set_contact_elim
subroutine, public hecmw_mat_set_contact_elim(hecMAT, contact_elim)
Definition: hecmw_matrix_misc.f90:461
m_fstr::myrank
integer(kind=kint) myrank
PARALLEL EXECUTION.
Definition: m_fstr.f90:96
m_solve_lineq_contact
This module provides functions to solve sparse system of \linear equitions in the case of contact ana...
Definition: hecmw_solver_contact.f90:9
m_solve_lineq_mumps_contact::solve_lineq_mumps_contact_init
subroutine, public solve_lineq_mumps_contact_init(hecMESH, hecMAT, hecLagMAT, is_sym)
Definition: solve_LINEQ_MUMPS_contact.f90:27
hecmw_matrix_misc::hecmw_mat_get_timelog
integer(kind=kint) function, public hecmw_mat_get_timelog(hecMAT)
Definition: hecmw_matrix_misc.f90:488
m_solve_lineq_direct_serial_lag
Definition: solve_LINEQ_direct_serial_lag.f90:5
m_solve_lineq_mkl_contact::solve_lineq_mkl_contact_init
subroutine, public solve_lineq_mkl_contact_init(hecMESH, is_sym)
Definition: solve_LINEQ_mkl.F90:30
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_matrix_misc::hecmw_mat_get_iterlog
integer(kind=kint) function, public hecmw_mat_get_iterlog(hecMAT)
Definition: hecmw_matrix_misc.f90:474
m_solve_lineq_direct_serial_lag::solve_lineq_serial_lag_hecmw_init
subroutine solve_lineq_serial_lag_hecmw_init(hecMAT, hecLagMAT, is_sym)
Definition: solve_LINEQ_direct_serial_lag.f90:17
m_solve_lineq_contact::solve_lineq_contact
subroutine, public solve_lineq_contact(hecMESH, hecMAT, hecLagMAT, conMAT, istat, rf, is_contact_active)
This subroutine.
Definition: hecmw_solver_contact.f90:74
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
m_solve_lineq_direct_serial_lag::solve_lineq_serial_lag_hecmw
subroutine solve_lineq_serial_lag_hecmw(hecMESH, hecMAT, hecLagMAT)
Definition: solve_LINEQ_direct_serial_lag.f90:28
m_solve_lineq_mkl_contact::solve_lineq_mkl_contact
subroutine, public solve_lineq_mkl_contact(hecMESH, hecMAT, hecLagMAT, istat, conMAT)
This subroutine executes the MKL solver.
Definition: solve_LINEQ_mkl.F90:56
m_solve_lineq_contact::solve_lineq_contact_init
subroutine, public solve_lineq_contact_init(hecMESH, hecMAT, hecLagMAT, is_sym)
This subroutine.
Definition: hecmw_solver_contact.f90:30
hecmw_matrix_misc
Definition: hecmw_matrix_misc.f90:6
m_solve_lineq_contact_elim::solve_lineq_contact_elim
subroutine, public solve_lineq_contact_elim(hecMESH, hecMAT, hecLagMAT, istat, conMAT, is_contact_active)
Definition: solve_LINEQ_contact_elim.f90:56
hecmw_matrix_misc::hecmw_mat_get_contact_elim
integer(kind=kint) function, public hecmw_mat_get_contact_elim(hecMAT)
Definition: hecmw_matrix_misc.f90:455
hecmw_util::hecmw_comm_get_size
integer(kind=kint) function hecmw_comm_get_size()
Definition: hecmw_util_f.F90:593
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
m_hecmw_comm_f::hecmw_update_r
subroutine hecmw_update_r(hecMESH, val, n, m)
Definition: hecmw_comm_f.F90:683
hecmw_util::hecmw_comm_get_comm
integer(kind=kint) function hecmw_comm_get_comm()
Definition: hecmw_util_f.F90:571
m_solve_lineq_contact_elim::solve_lineq_contact_elim_init
subroutine, public solve_lineq_contact_elim_init(hecMESH, hecMAT, hecLagMAT, is_sym)
Definition: solve_LINEQ_contact_elim.f90:34
hecmw_util::hecmwst_matrix
Definition: hecmw_util_f.F90:444