FrontISTR  5.7.1
Large-scale structural analysis program with finit element method
hecmw_solver_direct_ClusterMKL.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 !-------------------------------------------------------------------------------
6 
8  use hecmw_util
14 
15  private
17 
18  logical, save :: INITIALIZED = .false.
19  type (sparse_matrix), save :: spMAT
20 
21 contains
22 
23  subroutine hecmw_solve_direct_clustermkl(hecMESH,hecMAT)
24  implicit none
25  type (hecmwst_local_mesh), intent(in) :: hecmesh
26  type (hecmwst_matrix ), intent(inout) :: hecmat
27  integer(kind=kint) :: spmat_type
28  integer(kind=kint) :: spmat_symtype
29  integer(kind=kint) :: phase_start
30  integer(kind=kint) :: istat,myrank
31  real(kind=kreal) :: t1,t2
32 
33  t1=hecmw_wtime()
34 
35  call hecmw_mat_dump(hecmat, hecmesh)
36 
37  myrank=hecmw_comm_get_rank()
38 
39  if (initialized .and. hecmat%Iarray(98) .eq. 1) then
40  call sparse_matrix_finalize(spmat)
41  initialized = .false.
42  endif
43 
44  if (.not. initialized) then
45  spmat_type = sparse_matrix_type_csr
46  if (hecmat%symmetric) then
47  spmat_symtype = sparse_matrix_symtype_sym
48  else
49  spmat_symtype = sparse_matrix_symtype_asym
50  end if
51  call sparse_matrix_set_type(spmat, spmat_type, spmat_symtype)
52 
53  initialized = .true.
54  hecmat%Iarray(98) = 1
55  endif
56 
57  !* Flag to activate symbolic factorization: 1(yes) 0(no) hecMESH%Iarray(98)
58  !* Flag to activate numeric factorization: 1(yes) 0(no) hecMESH%Iarray(97)
59 
60  if (hecmat%Iarray(98) .eq. 1) then
61  ! ANALYSIS and FACTORIZATION
62  call sparse_matrix_hec_init_prof(spmat, hecmat, hecmesh)
63  call sparse_matrix_hec_set_vals(spmat, hecmat)
64  !call sparse_matrix_dump(spMAT)
65  phase_start = 1
66  hecmat%Iarray(98) = 0
67  hecmat%Iarray(97) = 0
68  endif
69  if (hecmat%Iarray(97) .eq. 1) then
70  ! FACTORIZATION
71  call sparse_matrix_hec_set_vals(spmat, hecmat)
72  !call sparse_matrix_dump(spMAT)
73  phase_start = 2
74  hecmat%Iarray(97) = 0
75  endif
76  call sparse_matrix_hec_set_rhs(spmat, hecmat)
77 
78  t2=hecmw_wtime()
79  if (myrank==0 .and. spmat%timelog > 0) &
80  write(*,'(A,f10.3)') ' [Cluster Pardiso]: Setup completed. time(sec)=',t2-t1
81 
82  ! SOLVE
83  call hecmw_clustermkl_wrapper(spmat, phase_start, hecmat%X, istat)
84 
85  call hecmw_mat_dump_solution(hecmat)
86 
87  end subroutine hecmw_solve_direct_clustermkl
88 
subroutine, public hecmw_mat_dump(hecMAT, hecMESH)
subroutine, public hecmw_mat_dump_solution(hecMAT)
This module provides linear equation solver interface for Cluster Pardiso.
subroutine, public hecmw_solve_direct_clustermkl(hecMESH, hecMAT)
I/O and Utility.
Definition: hecmw_util_f.F90:7
integer(kind=4), parameter kreal
integer(kind=kint) function hecmw_comm_get_rank()
real(kind=kreal) function hecmw_wtime()
This module provides linear equation solver interface for Cluster Pardiso.
subroutine, public hecmw_clustermkl_wrapper(spMAT, phase_start, solution, istat)
This module provides conversion routines between HEC data structure and DOF based sparse matrix struc...
subroutine, public sparse_matrix_hec_init_prof(spMAT, hecMAT, hecMESH)
subroutine, public sparse_matrix_hec_set_rhs(spMAT, hecMAT)
subroutine, public sparse_matrix_hec_set_vals(spMAT, hecMAT)
This module provides DOF based sparse matrix data structure (CSR and COO)
integer(kind=kint), parameter, public sparse_matrix_type_csr
integer(kind=kint), parameter, public sparse_matrix_symtype_sym
subroutine, public sparse_matrix_finalize(spMAT)
integer(kind=kint), parameter, public sparse_matrix_symtype_asym
subroutine, public sparse_matrix_set_type(spMAT, type, symtype)