FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hecmw_solver.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 contains
8 
9  subroutine hecmw_solve (hecMESH, hecMAT)
10 
11  use hecmw_util
20  use m_hecmw_comm_f
21  implicit none
22 
23  type (hecmwST_matrix), target :: hecMAT
24  type (hecmwST_local_mesh) :: hecMESH
25 
26  real(kind=kreal) :: resid
27  integer(kind=kint) :: i, myrank, NDOF
28  integer(kind=kint) :: imsg = 51
29  ndof=hecmat%NDOF
30 
31  !C ERROR CHECK
32  if(hecmw_solve_check_zerorhs(hecmesh, hecmat))then
33  hecmat%X = 0.0d0
34  return
35  endif
36 
37  select case(hecmat%Iarray(99))
38  !* Call Iterative Solver
39  case (1)
40  call hecmw_solve_iterative(hecmesh,hecmat)
41  !* Call Direct Solver
42  case(2:)
43  !C
44  !* Please note the following:
45  !* Flag to activate symbolic factorization: 1(yes) 0(no) hecMESH%Iarray(98)
46  !* Flag to activate numeric factorization: 1(yes) 0(no) hecMESH%Iarray(97)
47 
48  if (hecmat%Iarray(97) .gt. 1) hecmat%Iarray(97)=1
49 
50  call hecmw_mat_set_flag_converged(hecmat, 0)
51  call hecmw_mat_set_flag_diverged(hecmat, 0)
52 
53  if (hecmat%Iarray(2) .eq. 102) then
54  if(hecmesh%PETOT.GT.1) then
55  call hecmw_solve_direct_clustermkl(hecmesh, hecmat)
56  else
57  call hecmw_solve_direct_mkl(hecmesh,hecmat)
58  endif
59  elseif (hecmat%Iarray(2) .eq. 104) then
60  call hecmw_solve_direct_mumps(hecmesh, hecmat)
61  else
62  if(hecmesh%PETOT.GT.1) then
63  call hecmw_solve_direct_parallel(hecmesh,hecmat,imsg)
64  else
65  call hecmw_solve_direct(hecmesh,hecmat,imsg)
66  endif
67  endif
68 
69  resid=hecmw_rel_resid_l2(hecmesh,hecmat)
70  myrank=hecmw_comm_get_rank()
71  if (myrank==0) then
72  if (hecmat%Iarray(21) > 0 .or. hecmat%Iarray(22) > 0) then
73  write(*,"(a,1pe12.5)")'### Relative residual =', resid
74  endif
75  if( resid >= 1.0d-8) then
76  write(*,"(a)")'### Relative residual exceeded 1.0d-8---Direct Solver### '
77  ! stop
78  endif
79  endif
80  if (resid < hecmw_mat_get_resid(hecmat)) then
81  call hecmw_mat_set_flag_converged(hecmat, 1)
82  endif
83  !C
84  end select
85 
86  end subroutine hecmw_solve
87 
88  subroutine hecmw_substitute_solver(hecMESH, hecMATorig, NDOF)
89 
90  use hecmw_util
96  type (hecmwST_local_mesh) :: hecMESH
97  type (hecmwST_matrix) :: hecMATorig
98  type (hecmwST_matrix),pointer :: hecMAT => null()
99  integer(kind=kint) NDOF
100  if (ndof == hecmatorig%NDOF) then
101  call hecmw_clone_matrix(hecmatorig,hecmat)
102  else if (ndof < hecmatorig%NDOF) then
104  end if
105  ! select case(NDOF)
106  ! case(1)
107  ! call hecmw_solve_11(hecMESH,hecMAT)
108  ! case(2)
109  ! call hecmw_solve_22(hecMESH,hecMAT)
110  ! case(3)
111  ! call hecmw_solve_33(hecMESH,hecMAT)
112  ! case(4)
113  ! call hecmw_solve_iterative(hecMESH,hecMAT)
114  ! case(5)
115  ! call hecmw_solve_iterative(hecMESH,hecMAT)
116  ! case(6)
117  ! call hecmw_solve_66(hecMESH,hecMAT)
118  ! case(7:)
119  ! call hecmw_solve_iterative(hecMESH,hecMAT)
120  ! end select
121  !call hecmw_solve_direct_MUMPS(hecMESH, hecMAT)
122  call hecmw_solve_iterative(hecmesh,hecmat)
123  if (ndof /= hecmatorig%NDOF) then
124  call hecmw_vector_contract(hecmatorig,hecmat,ndof)
125  end if
126  end subroutine hecmw_substitute_solver
127 
128 end module hecmw_solver
hecmw_solver::hecmw_solve
subroutine hecmw_solve(hecMESH, hecMAT)
Definition: hecmw_solver.f90:10
hecmw_solver_direct_mkl
This module provides linear equation solver interface for Pardiso.
Definition: hecmw_solver_direct_MKL.f90:7
hecmw_util::hecmw_abort
subroutine hecmw_abort(comm)
Definition: hecmw_util_f.F90:534
hecmw_util::hecmw_clone_matrix
subroutine hecmw_clone_matrix(hecMATorig, hecMAT)
Definition: hecmw_util_f.F90:916
hecmw_solver_iterative
Definition: hecmw_solver_Iterative.f90:6
hecmw_solver_las
Definition: hecmw_solver_las.f90:6
hecmw_solver_iterative::hecmw_solve_iterative
subroutine hecmw_solve_iterative(hecMESH, hecMAT)
Definition: hecmw_solver_Iterative.f90:14
hecmw_matrix_misc::hecmw_mat_get_resid
real(kind=kreal) function, public hecmw_mat_get_resid(hecMAT)
Definition: hecmw_matrix_misc.f90:672
hecmw_solver_direct_mkl::hecmw_solve_direct_mkl
subroutine, public hecmw_solve_direct_mkl(hecMESH, hecMAT)
Definition: hecmw_solver_direct_MKL.f90:24
hecmw_solver_direct_parallel::hecmw_solve_direct_parallel
subroutine, public hecmw_solve_direct_parallel(hecMESH, hecMAT, ii)
Definition: hecmw_solver_direct_parallel.F90:82
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
hecmw_matrix_misc::hecmw_mat_set_flag_diverged
subroutine, public hecmw_mat_set_flag_diverged(hecMAT, flag_diverged)
Definition: hecmw_matrix_misc.f90:625
m_hecmw_comm_f
Definition: hecmw_comm_f.F90:6
hecmw_solver_direct_mumps
This module provides linear equation solver interface for MUMPS.
Definition: hecmw_solver_direct_MUMPS.f90:6
hecmw_solver_direct_parallel
Definition: hecmw_solver_direct_parallel.F90:6
hecmw_solver_direct::hecmw_solve_direct
subroutine, public hecmw_solve_direct(hecMESH, hecMAT, Ifmsg)
HECMW_SOLVE_DIRECT is a program for the matrix solver.
Definition: hecmw_solver_direct.f90:61
hecmw_solver::hecmw_substitute_solver
subroutine hecmw_substitute_solver(hecMESH, hecMATorig, NDOF)
Definition: hecmw_solver.f90:89
hecmw_solver_las::hecmw_rel_resid_l2
real(kind=kreal) function, public hecmw_rel_resid_l2(hecMESH, hecMAT, COMMtime)
Definition: hecmw_solver_las.f90:125
hecmw_solver_direct
HECMW_SOLVE_DIRECT is a program for the matrix direct solver.
Definition: hecmw_solver_direct.f90:10
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
hecmw_matrix_misc
Definition: hecmw_matrix_misc.f90:6
hecmw_solver_iterative::hecmw_solve_check_zerorhs
logical function hecmw_solve_check_zerorhs(hecMESH, hecMAT)
Definition: hecmw_solver_Iterative.f90:282
hecmw_solver_direct_mumps::hecmw_solve_direct_mumps
subroutine, public hecmw_solve_direct_mumps(hecMESH, hecMAT)
Definition: hecmw_solver_direct_MUMPS.f90:23
hecmw_matrix_misc::hecmw_mat_set_flag_converged
subroutine, public hecmw_mat_set_flag_converged(hecMAT, flag_converged)
Definition: hecmw_matrix_misc.f90:613
hecmw_util::hecmw_comm_get_rank
integer(kind=kint) function hecmw_comm_get_rank()
Definition: hecmw_util_f.F90:582
hecmw_solver_direct_clustermkl::hecmw_solve_direct_clustermkl
subroutine, public hecmw_solve_direct_clustermkl(hecMESH, hecMAT)
Definition: hecmw_solver_direct_ClusterMKL.f90:24
hecmw_util::hecmw_comm_get_comm
integer(kind=kint) function hecmw_comm_get_comm()
Definition: hecmw_util_f.F90:571
hecmw_util::hecmw_vector_contract
subroutine hecmw_vector_contract(hecMATorig, hecMAT, NDOF)
Definition: hecmw_util_f.F90:1048
hecmw_solver
Definition: hecmw_solver.f90:6
hecmw_solver_direct_clustermkl
This module provides linear equation solver interface for Cluster Pardiso.
Definition: hecmw_solver_direct_ClusterMKL.f90:7