SPEED
MAKE_SYSTEM_POSITION_LGLNODES.f90
Go to the documentation of this file.
1! Copyright (C) 2012 The SPEED FOUNDATION
2! Author: Ilario Mazzieri
3!
4! This file is part of SPEED.
5!
6! SPEED is free software; you can redistribute it and/or modify it
7! under the terms of the GNU Affero General Public License as
8! published by the Free Software Foundation, either version 3 of the
9! License, or (at your option) any later version.
10!
11! SPEED is distributed in the hope that it will be useful, but
12! WITHOUT ANY WARRANTY; without even the implied warranty of
13! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14! Affero General Public License for more details.
15!
16! You should have received a copy of the GNU Affero General Public License
17! along with SPEED. If not, see <http://www.gnu.org/licenses/>.
18
26
27subroutine make_system_position_lglnodes(nn_loc , local_node_num, cs_nnz_loc, cs_loc, &
28 xs_loc, ys_loc, zs_loc, nsdof, sys_label, &
29 x_system_lst, y_system_lst, z_system_lst, &
30 mpi_id, mpi_comm, mpi_np)
31
33
34 implicit none
35
36 include 'SPEED.MPI'
37
38 integer*4 :: nn_loc, nsdof, cs_nnz_loc
39 integer*4 :: mpi_np, mpi_id, mpi_comm, mpi_ierr
40 integer*4 :: i, j, max_nsdof_per_node
41 integer*4, dimension(0:cs_nnz_loc) :: cs_loc
42 integer*4, dimension(nn_loc) :: node_counter_loc
43 integer*4, dimension(nsdof) :: node_counter_sdof_dum
44 integer*4, dimension(nsdof) :: node_sys_loc, sys_label, nearestnode_globnum, nearestnode_locnum
45 integer*4, dimension(nsdof*mpi_np) :: node_sys_glo
46 integer*4, dimension(nn_loc) :: local_node_num
47
48 !integer*4, dimension(:,:), allocatable :: locnode_buildID_map !INTENT(OUT)
49
50 real*8, dimension(nn_loc) :: xs_loc, ys_loc, zs_loc
51 real*8, dimension(nsdof) :: x_system_lst, y_system_lst, z_system_lst
52 real*8, dimension(nsdof) :: dist_system_loc
53 real*8, dimension(nsdof*mpi_np) :: dist_system_glo
54
55 node_counter_loc = 0
56
57 ! Find if LGL_gob_num is inside the local partition
58 ! Finding Nearest Node to SDOF Oscillator across all partitions - in Global Numeration
59 do i = 1,nsdof
60 call get_nearest_node(nn_loc, xs_loc, ys_loc, zs_loc, x_system_lst(i), y_system_lst(i), z_system_lst(i), &
61 node_sys_loc(i), dist_system_loc(i))
62 node_sys_loc(i) = local_node_num(node_sys_loc(i))
63 enddo
64
65 call mpi_barrier(mpi_comm, mpi_ierr)
66
67 call mpi_allgather(dist_system_loc, nsdof, speed_double, dist_system_glo, nsdof, speed_double, mpi_comm, mpi_ierr)
68 call mpi_allgather(node_sys_loc, nsdof, speed_integer, node_sys_glo, nsdof, speed_integer, mpi_comm, mpi_ierr)
69
70 node_sys_loc = 0
71
72 call get_minvalues(node_sys_glo, dist_system_glo, nsdof*mpi_np, node_sys_loc, nsdof, mpi_np)
73 call mpi_barrier(mpi_comm, mpi_ierr)
74
75 ! Finding if this nearest node is common in multiple partitions
76 node_counter_sdof_dum = 0
77 do i = 1, nsdof
78 nearestnode_globnum(i) = node_sys_glo(node_sys_loc(i))
79
80 call get_indloc_from_indglo(local_node_num, nn_loc, nearestnode_globnum(i), nearestnode_locnum(i))
81
82 if (nearestnode_locnum(i).ne.0) then
83 node_counter_sdof_dum(i) = 1;
84 node_counter_loc(nearestnode_locnum(i)) = node_counter_loc(nearestnode_locnum(i)) + 1
85 endif
86 enddo
87
88 ! Finiding LGL nodes that are shared by different partitions
89 allocate(node_counter_sdof(nsdof))
90 call mpi_barrier(mpi_comm, mpi_ierr)
91 call mpi_allreduce(node_counter_sdof_dum, node_counter_sdof, nsdof, speed_integer, mpi_sum, mpi_comm, mpi_ierr)
92
93 max_nsdof_per_node = maxval(node_counter_loc) ! num on SDOFs allocated to a LGL node, id of 1st SDOF, id of 2nd SDOF etc...
94 allocate(locnode_buildid_map(nn_loc, max_nsdof_per_node + 1))
96
97 ! Mapping SDOFs to nearest LGL nodes;
98 do i = 1, nsdof
99 nearestnode_globnum(i) = node_sys_glo(node_sys_loc(i))
100
101 call get_indloc_from_indglo(local_node_num, nn_loc, nearestnode_globnum(i), nearestnode_locnum(i))
102
103 if (nearestnode_locnum(i).ne.0) then
104 locnode_buildid_map(nearestnode_locnum(i), 1) = locnode_buildid_map(nearestnode_locnum(i), 1) + 1
105 locnode_buildid_map(nearestnode_locnum(i), locnode_buildid_map(nearestnode_locnum(i), 1)+1 ) = sys_label(i)
106 ! write(*,*) mpi_id, i, sys_label(i),'; node_counter = ', node_counter_sdof, '; coords = ', xs_loc(nearestnode_locnum(i)), &
107 ! ys_loc(nearestnode_locnum(i)), zs_loc(nearestnode_locnum(i))
108 endif
109 enddo
110
subroutine get_indloc_from_indglo(local_el, nel_loc, ie, ic)
Returns local id from global id.
subroutine get_minvalues(i_glo, v_glo, n_glo, i_loc, n_loc, np)
Computes positions of minimum values.
subroutine get_nearest_node(n, xs, ys, zs, xt, yt, zt, nt, dist_min)
Computes the nearest node with respet to (xt,yt,zt).
subroutine make_system_position_lglnodes(nn_loc, local_node_num, cs_nnz_loc, cs_loc, xs_loc, ys_loc, zs_loc, nsdof, sys_label, x_system_lst, y_system_lst, z_system_lst, mpi_id, mpi_comm, mpi_np)
Reads oscillator position from SYS.input, and finds to which local-LGL node it should be applied.
Contains SPEED PARAMETERS used in (SPEED, READ_INPUT_FILES, MAKE_PARTION_AND_MPI_FILES,...
Definition MODULES.f90:207
integer *4, dimension(:,:), allocatable locnode_buildid_map
Definition MODULES.f90:353
real *8, dimension(:), allocatable y_system_lst
Definition MODULES.f90:417
real *8, dimension(:), allocatable dist_system_glo
Definition MODULES.f90:417
real *8, dimension(:), allocatable x_system_lst
Definition MODULES.f90:417
real *8, dimension(:), allocatable z_system_lst
Definition MODULES.f90:417
integer *4, dimension(:), allocatable node_counter_sdof
Definition MODULES.f90:355