SPEED
SETUP_DG4MPI.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine setup_dg4mpi (el_new, nel_dg_loc, nel_dg_glo, mpi_id, mpi_comm, mpi_np, cs_nnz_loc, cs_loc, ne_loc, local_el_num, local_n_num, nn_loc, total_cs_nnz_mpi, total_el, mpi_file)
 Writes DGCSXXXXXX.mpi files for DG connectivity.
 

Function/Subroutine Documentation

◆ setup_dg4mpi()

subroutine setup_dg4mpi ( type(el4loop), dimension(nel_dg_loc), intent(in)  el_new,
integer*4  nel_dg_loc,
integer*4  nel_dg_glo,
integer*4  mpi_id,
integer*4  mpi_comm,
integer*4  mpi_np,
integer*4  cs_nnz_loc,
integer*4, dimension(0:cs_nnz_loc)  cs_loc,
integer*4  ne_loc,
integer*4, dimension(ne_loc)  local_el_num,
integer*4, dimension(nn_loc)  local_n_num,
integer*4  nn_loc,
integer*4, intent(out)  total_cs_nnz_mpi,
integer*4, intent(out)  total_el,
character*70  mpi_file 
)

Writes DGCSXXXXXX.mpi files for DG connectivity.

Author
Ilario Mazzieri
Date
September, 2013
Version
1.0
Parameters
[in]el_newstructure for DG faces
[in]nel_dg_locnumber of local dg elements
[in]nel_dg_glonumber of global dg elements
[in]mpi_idmpi processor identity
[in]mpi_commmpi communicator
[in]mpi_npnumber of mpi processes
[in]cs_nnz_loclength of cs_loc
[in]cs_loclocal spectral connectivity vector
[in]nn_locnumber of local nodes
[in]local_n_numlocal node numeration
[in]local_el_numlocal element numeration
[in]mpi_filefolder where to write *.mpi file
[out]DGCSXXXXXX.mpiDG connectivity vector for MPI routines (XXXXXX -> processor number)
[out]total_cs_nnz_mpitotal number of DG nodes for mpi routines
[out]total_eltotal number of DG elements for mpi routines

Definition at line 39 of file SETUP_DG4MPI.f90.

43
44 use max_var
45 use dgjump
46
47 implicit none
48
49
50 include 'SPEED.MPI'
51
52 type(el4loop), dimension(nel_dg_loc), intent(in) :: el_new
53
54 character*70 :: file_mpi, mpi_file, file_mpi_new
55 character*70 :: filename
56
57 integer*4 :: nel_dg_loc, cs_nnz_loc, ne_loc, nel_dg_glo, mpi_id, mpi_comm, mpi_np, nn_loc
58 integer*4 :: mpierror, cs_nnz_mpi, nofel_mpi, el_dg
59 integer*4 :: i, j, k, ie, max_dime, dime_mpi, unit_mpi
60 integer*4, intent(out) :: total_cs_nnz_mpi, total_el
61
62 integer*4, dimension(:), allocatable :: cs_mpi
63 integer*4, dimension(0:cs_nnz_loc) :: cs_loc
64 integer*4, dimension(ne_loc) :: local_el_num
65 integer*4, dimension(nn_loc) :: local_n_num
66
67 cs_nnz_mpi = 0
68 nofel_mpi = 1
69
70 do i = 1, nel_dg_loc
71 if (i .eq. 1) then
72 cs_nnz_mpi = cs_nnz_mpi + el_new(i)%deg**3 + 1 + 1
73 nofel_mpi = nofel_mpi + 1
74
75 elseif(el_new(i)%ind .ne. el_new(i-1)%ind) then
76 cs_nnz_mpi = cs_nnz_mpi + el_new(i)%deg**3 + 1 + 1
77 nofel_mpi = nofel_mpi + 1
78 endif
79 enddo
80
81 allocate(cs_mpi(0:cs_nnz_mpi))
82 cs_mpi = 0
83 cs_mpi(0) = nofel_mpi
84 k = 0
85
86
87 do i = 1, nel_dg_loc
88
89 if (i .eq. 1) then
90
91
92 call get_indloc_from_indglo(local_el_num, ne_loc, el_new(i)%ind, ie)
93 k = 1
94 cs_mpi(k) = cs_mpi(k-1) + el_new(i)%deg**3 + 1
95 cs_mpi(cs_mpi(k-1)) = el_new(i)%ind
96
97 do j= 1, el_new(i)%deg**3
98 cs_mpi(cs_mpi(k-1) + j) = local_n_num(cs_loc(cs_loc(ie-1) + j))
99 enddo
100
101 elseif(el_new(i)%ind .ne. el_new(i-1)%ind) then
102
103 call get_indloc_from_indglo(local_el_num, ne_loc, el_new(i)%ind, ie)
104 k = k + 1
105
106 cs_mpi(k) = cs_mpi(k-1) + el_new(i)%deg**3 + 1
107 cs_mpi(cs_mpi(k-1)) = el_new(i)%ind
108
109 do j= 1, el_new(i)%deg**3
110 cs_mpi(cs_mpi(k-1) + j) = local_n_num(cs_loc(cs_loc(ie-1) + j))
111 enddo
112
113 endif
114 enddo
115
116
117! write(*,*) 'dopo', cs_mpi
118
119
120 file_mpi = 'DGCS000000.mpi'
121 unit_mpi = 40
122
123 if (mpi_id .lt. 10) then
124 write(file_mpi(10:10),'(i1)') mpi_id
125 else if (mpi_id .lt. 100) then
126 write(file_mpi(9:10),'(i2)') mpi_id
127 else if (mpi_id .lt. 1000) then
128 write(file_mpi(8:10),'(i3)') mpi_id
129 else if (mpi_id .lt. 10000) then
130 write(file_mpi(7:10),'(i4)') mpi_id
131 else if (mpi_id .lt. 100000) then
132 write(file_mpi(6:10),'(i5)') mpi_id
133 else if (mpi_id .lt. 1000000) then
134 write(file_mpi(5:10),'(i6)') mpi_id
135 endif
136
137 if(len_trim(mpi_file) .ne. 70) then
138 file_mpi_new = mpi_file(1:len_trim(mpi_file)) // '/' // file_mpi
139 else
140 file_mpi_new = file_mpi
141 endif
142
143 open(unit_mpi,file=file_mpi_new)
144 write(unit_mpi,*) cs_nnz_mpi
145 do i = 0, cs_nnz_mpi
146 write(unit_mpi,*) cs_mpi(i)
147 enddo
148 close(unit_mpi)
149
150
151 deallocate(cs_mpi)
152
153 call mpi_barrier(mpi_comm, mpierror)
154
155 total_cs_nnz_mpi = 0
156 total_el = 0
157
158 do i = 1, mpi_np
159 file_mpi = 'DGCS000000.mpi'
160 unit_mpi = 40
161
162 if (i-1 .lt. 10) then
163 write(file_mpi(10:10),'(i1)') i-1
164 else if (i-1 .lt. 100) then
165 write(file_mpi(9:10),'(i2)') i-1
166 else if (i-1 .lt. 1000) then
167 write(file_mpi(8:10),'(i3)') i-1
168 else if (i-1 .lt. 10000) then
169 write(file_mpi(7:10),'(i4)') i-1
170 else if (i-1 .lt. 100000) then
171 write(file_mpi(6:10),'(i5)') i-1
172 else if (i-1 .lt. 1000000) then
173 write(file_mpi(5:10),'(i6)') i-1
174 endif
175
176 if(len_trim(mpi_file) .ne. 70) then
177 file_mpi_new = mpi_file(1:len_trim(mpi_file)) // '/' // file_mpi
178 else
179 file_mpi_new = file_mpi
180 endif
181
182 open(unit_mpi,file=file_mpi_new)
183 read(unit_mpi,*) cs_nnz_mpi
184 read(unit_mpi,*) el_dg
185
186 total_cs_nnz_mpi = total_cs_nnz_mpi + cs_nnz_mpi
187 total_el = total_el + el_dg -1
188 close(unit_mpi)
189
190 enddo
191
192 total_el = total_el + 1
193
194 return
195
196
subroutine get_indloc_from_indglo(local_el, nel_loc, ie, ic)
Returns local id from global id.
Contains structure for jump matrices.
Definition MODULES.f90:155
Set maximal bounds.
Definition MODULES.f90:54

References get_indloc_from_indglo().

Here is the call graph for this function: