SFEMaNS  version 4.1 (work in progress)
Reference documentation for SFEMaNS
 All Classes Files Functions Variables Groups Pages
prep_mesh_periodic.f90
Go to the documentation of this file.
1 !
2 !Authors Jean-Luc Guermond, Copyrights 1996
3 !
4 MODULE periodic
5 
6 
7  IMPLICIT NONE
8 
9  PUBLIC :: prep_periodic, &
14  PRIVATE
15 
16 CONTAINS
17 
18  !jan 29 2007
19  SUBROUTINE prep_periodic(my_periodic, mesh, periodic)
20  !=========================================
21  USE def_type_mesh
23  IMPLICIT NONE
24  TYPE(periodic_data), INTENT(IN) :: my_periodic
25  TYPE(mesh_type) :: mesh
26  TYPE(periodic_type) :: periodic
27  INTEGER, DIMENSION(:), POINTER :: list_loc, perlist_loc, list_dom, perlist_dom
28  INTEGER :: n, side1, side2, n_b, nx, i
29  REAL(KIND=8), DIMENSION(:), POINTER :: e
30 
31  WRITE (*,*) 'Loading periodic-data file ...'
32 
33  IF (mesh%np == 0) THEN
34  WRITE(*,*) 'no mesh on this proc'
35  RETURN
36  END IF
37 
38  ALLOCATE(e(SIZE(my_periodic%vect_e,1)))
39 
40  periodic%n_bord = my_periodic%nb_periodic_pairs
41  IF (periodic%n_bord .GT. 20) THEN
42  WRITE(*,*) 'PREP_MESH_PERIODIC: trop de bords periodiques'
43  stop
44  END IF
45 
46  DO n= 1, periodic%n_bord
47 
48  side1 = my_periodic%list_periodic(1,n)
49  side2 = my_periodic%list_periodic(2,n)
50  e = my_periodic%vect_e(:,n)
51 
52  CALL list_periodic(mesh%np, mesh%jjs, mesh%sides, mesh%rr, side1, side2, e, &
53  list_loc, perlist_loc)
54 
55  n_b = SIZE(list_loc)
56  ALLOCATE(list_dom(n_b), perlist_dom(n_b))
57  nx = 0
58  DO i = 1, n_b
59  IF (max(list_loc(i),perlist_loc(i)) .LE. mesh%dom_np) THEN
60  nx = nx+1
61  list_dom(nx) = list_loc(i)
62  perlist_dom(nx) = perlist_loc(i)
63  ELSE IF (min(list_loc(i),perlist_loc(i)) .LE. mesh%dom_np) THEN
64  WRITE(*,*) 'BUG in prep_periodic'
65  stop
66  END IF
67  END DO
68  IF (n_b /= nx) WRITE(*,*) 'WARNING, I have removed', n_b-nx, ' periodic pairs in prep_periodic'
69  n_b = nx
70 
71  ALLOCATE (periodic%list(n)%DIL(n_b), periodic%perlist(n)%DIL(n_b))
72  periodic%list(n)%DIL = list_dom(1:n_b)
73  periodic%perlist(n)%DIL = perlist_dom(1:n_b)
74 
75  DEALLOCATE(list_loc,perlist_loc, list_dom, perlist_dom)
76  END DO
77 
78  DEALLOCATE(e)
79 
80  WRITE (*,*) 'Treatment of periodic-data done'
81 
82  END SUBROUTINE prep_periodic
83 
84  SUBROUTINE prep_periodic_bloc(my_periodic, mesh, periodic, nb_bloc)
85  !=========================================
87  USE def_type_mesh
88  IMPLICIT NONE
89  TYPE(periodic_data), INTENT(IN) :: my_periodic
90  TYPE(mesh_type) :: mesh
91  TYPE(periodic_type) :: periodic
92  INTEGER, INTENT(IN) :: nb_bloc
93  INTEGER, DIMENSION(:), POINTER :: list_loc, perlist_loc, list_dom, perlist_dom
94  INTEGER :: n, side1, side2, nsize, n_b
95  INTEGER :: k, k_deb, k_fin, nx, i
96  REAL(KIND=8), DIMENSION(2) :: e
97 
98  WRITE (*,*) 'Loading periodic-data file ...'
99 
100  IF (mesh%np == 0) THEN
101  WRITE(*,*) 'no mesh on this proc'
102  RETURN
103  END IF
104 
105  periodic%n_bord = my_periodic%nb_periodic_pairs
106  IF (periodic%n_bord .GT. 20) THEN
107  WRITE(*,*) 'PREP_MESH_PERIODIC: trop de bords periodiques'
108  stop
109  END IF
110 
111  DO n= 1, periodic%n_bord
112 
113  side1 = my_periodic%list_periodic(1,n)
114  side2 = my_periodic%list_periodic(2,n)
115  e = my_periodic%vect_e(:,n)
116 
117  CALL list_periodic(mesh%np, mesh%jjs, mesh%sides, mesh%rr, side1, side2, e, &
118  list_loc, perlist_loc)
119 
120  !n_b = SIZE(list_loc)
121  n_b = SIZE(perlist_loc)
122  ALLOCATE(list_dom(n_b), perlist_dom(n_b))
123  nx = 0
124  DO i = 1, n_b
125  IF (max(list_loc(i),perlist_loc(i)) .LE. mesh%dom_np) THEN
126  nx = nx+1
127  list_dom(nx) = list_loc(i)
128  perlist_dom(nx) = perlist_loc(i)
129  ELSE IF (min(list_loc(i),perlist_loc(i)) .LE. mesh%dom_np) THEN
130  WRITE(*,*) 'BUG in prep_periodic_bloc'
131  stop
132  END IF
133  END DO
134  IF (n_b /= nx) WRITE(*,*) 'WARNING, I have removed', n_b-nx, ' periodic pairs in prep_periodic_bloc'
135  n_b = nx
136 
137  nsize = nb_bloc*n_b !SIZE(list_loc) !n_b
138 
139  ALLOCATE(periodic%list(n)%DIL(nsize), periodic%perlist(n)%DIL(nsize))
140 
141  DO k = 1, nb_bloc
142  k_deb=(k-1)*n_b+1
143  k_fin=k*n_b
144  periodic%list(n)%DIL(k_deb:k_fin) = list_dom(1:n_b) + (k-1)*mesh%dom_np ! First bloc
145  periodic%perlist(n)%DIL(k_deb:k_fin) = perlist_dom(1:n_b) + (k-1)*mesh%dom_np ! First bloc
146  END DO
147 
148  DEALLOCATE(list_loc,perlist_loc, list_dom, perlist_dom)
149 
150  END DO
151 
152  WRITE (*,*) 'Treatment of periodic-data done'
153 
154  END SUBROUTINE prep_periodic_bloc
155  !jan 29 2007
156 
157  !JLG+FL/Feb 2 2010
158  SUBROUTINE prep_periodic_h_p_phi_bc(my_periodic, H_mesh, pmag_mesh, phi_mesh, H_p_phi_per)
159  USE chaine_caractere
160  USE def_type_mesh
161  IMPLICIT NONE
162  TYPE(periodic_data), INTENT(IN) :: my_periodic
163  TYPE(mesh_type) :: h_mesh, pmag_mesh, phi_mesh
164  TYPE(periodic_type) :: h_p_phi_per
165  INTEGER, DIMENSION(:), POINTER :: b_list_loc, b_perlist_loc, &
166  e_list_loc, e_perlist_loc, p_list_loc, p_perlist_loc
167  INTEGER, DIMENSION(:), POINTER :: b_list_dom, b_perlist_dom, &
168  e_list_dom, e_perlist_dom, p_list_dom, p_perlist_dom
169  INTEGER :: n, side1, side2, nsize, n_b, n_e, n_p, nx, i
170  REAL(KIND=8), DIMENSION(2) :: e
171 
172  WRITE (*,*) 'Loading periodic-data file ...'
173 
174  h_p_phi_per%n_bord = my_periodic%nb_periodic_pairs
175  IF (h_p_phi_per%n_bord .GT. 20) THEN
176  WRITE(*,*) 'PREP_MESH_PERIODIC: trop de bords periodiques'
177  stop
178  END IF
179 
180  DO n= 1, h_p_phi_per%n_bord
181 
182  side1 = my_periodic%list_periodic(1,n)
183  side2 = my_periodic%list_periodic(2,n)
184  e = my_periodic%vect_e(:,n)
185 
186 
187  IF (h_mesh%np == 0) THEN
188  n_b = 0
189  ELSE
190  CALL list_periodic(h_mesh%np, h_mesh%jjs, h_mesh%sides, h_mesh%rr, side1, side2, e, &
191  b_list_loc, b_perlist_loc)
192  n_b = SIZE(b_list_loc)
193  ALLOCATE(b_list_dom(n_b), b_perlist_dom(n_b))
194  nx = 0
195  DO i = 1, n_b
196  IF (max(b_list_loc(i),b_perlist_loc(i)) .LE. h_mesh%dom_np) THEN
197  nx = nx+1
198  b_list_dom(nx) = b_list_loc(i)
199  b_perlist_dom(nx) = b_perlist_loc(i)
200  ELSE IF (min(b_list_loc(i),b_perlist_loc(i)) .LE. h_mesh%dom_np) THEN
201  WRITE(*,*) 'BUG in prep_periodic_H_p_phi_bc (H)'
202  stop
203  END IF
204  END DO
205  IF (n_b /= nx) WRITE(*,*) 'WARNING, I have removed', n_b-nx, ' periodic pairs on H'
206  n_b = nx
207  END IF
208 
209  IF (pmag_mesh%np == 0) THEN
210  n_p = 0
211  ELSE
212  CALL list_periodic(pmag_mesh%np, pmag_mesh%jjs, pmag_mesh%sides, pmag_mesh%rr, side1, side2, e, &
213  p_list_loc, p_perlist_loc)
214  n_p = SIZE(p_list_loc)
215  ALLOCATE(p_list_dom(n_p), p_perlist_dom(n_p))
216  nx = 0
217  DO i = 1, n_p
218  IF (max(p_list_loc(i),p_perlist_loc(i)) .LE. pmag_mesh%dom_np) THEN
219  nx = nx+1
220  p_list_dom(nx) = p_list_loc(i)
221  p_perlist_dom(nx) = p_perlist_loc(i)
222  ELSE IF (min(p_list_loc(i),p_perlist_loc(i)) .LE. pmag_mesh%dom_np) THEN
223  WRITE(*,*) 'BUG in prep_periodic_H_p_phi_bc (pmag) '
224  stop
225  END IF
226  END DO
227  IF (n_p /= nx) WRITE(*,*) 'WARNING, I have removed', n_p-nx, ' periodic pairs on pmag'
228  n_p = nx
229  END IF
230 
231  IF (phi_mesh%np==0) THEN
232  n_e = 0
233  ELSE
234  CALL list_periodic(phi_mesh%np, phi_mesh%jjs, phi_mesh%sides, phi_mesh%rr, side1, side2, e, &
235  e_list_loc, e_perlist_loc)
236  n_e = SIZE(e_list_loc)
237  ALLOCATE(e_list_dom(n_e), e_perlist_dom(n_e))
238  nx = 0
239  DO i = 1, n_e
240  IF (max(e_list_loc(i),e_perlist_loc(i)) .LE. phi_mesh%dom_np) THEN
241  nx = nx+1
242  e_list_dom(nx) = e_list_loc(i)
243  e_perlist_dom(nx) = e_perlist_loc(i)
244  ELSE IF (min(e_list_loc(i),e_perlist_loc(i)) .LE. phi_mesh%dom_np) THEN
245  WRITE(*,*) 'BUG in prep_periodic_H_p_phi_bc (phi) '
246  stop
247  END IF
248  END DO
249  IF (n_e /= nx) WRITE(*,*) 'WARNING, I have removed', n_e-nx, ' periodic pairs on phi'
250  n_e = nx
251  END IF
252 
253  !n_b = SIZE(b_list_loc)
254  !n_p = SIZE(p_list_loc)
255  !n_e = SIZE(e_list_loc)
256  nsize = 3*n_b + n_p + n_e
257 
258  ALLOCATE(h_p_phi_per%list(n)%DIL(nsize), h_p_phi_per%perlist(n)%DIL(nsize))
259 
260  IF (n_b /=0) THEN
261 
262  h_p_phi_per%list(n)%DIL(1:n_b) = b_list_dom(1:n_b) ! First block
263  h_p_phi_per%perlist(n)%DIL(1:n_b) = b_perlist_dom(1:n_b) ! First block
264 
265  h_p_phi_per%list(n)%DIL(n_b+1:2*n_b) = b_list_dom(1:n_b) + h_mesh%dom_np ! Second block
266  h_p_phi_per%perlist(n)%DIL(n_b+1:2*n_b) = b_perlist_dom(1:n_b) + h_mesh%dom_np ! Second block
267 
268  h_p_phi_per%list(n)%DIL(2*n_b+1:3*n_b) = b_list_dom(1:n_b) + 2*h_mesh%dom_np ! Third block
269  h_p_phi_per%perlist(n)%DIL(2*n_b+1:3*n_b) = b_perlist_dom(1:n_b) + 2*h_mesh%dom_np ! Third block
270 
271  END IF
272 
273  IF (n_p /=0) THEN
274  h_p_phi_per%list(n)%DIL(3*n_b+1:3*n_b+n_p) = p_list_dom(1:n_p) + 3*h_mesh%dom_np ! Forth block
275  h_p_phi_per%perlist(n)%DIL(3*n_b+1:3*n_b+n_p) = p_perlist_dom(1:n_p) + 3*h_mesh%dom_np ! Forth block
276  END IF
277 
278  IF (n_e/=0) THEN
279  h_p_phi_per%list(n)%DIL(3*n_b+n_p+1:) = e_list_dom(1:n_e) + 3*h_mesh%dom_np + pmag_mesh%dom_np ! Fourth block
280  h_p_phi_per%perlist(n)%DIL(3*n_b+n_p+1:) = e_perlist_dom(1:n_e) + 3*h_mesh%dom_np + pmag_mesh%dom_np ! Fourth block
281  END IF
282 
283  IF (ASSOCIATED(b_list_loc)) nullify(b_list_loc, b_perlist_loc, b_list_dom, b_perlist_dom)
284  IF (ASSOCIATED(p_list_loc)) nullify(p_list_loc, p_perlist_loc, p_list_dom, p_perlist_dom)
285  IF (ASSOCIATED(e_list_loc)) nullify(e_list_loc, e_perlist_loc, e_list_dom, e_perlist_dom)
286 
287  END DO
288 
289  WRITE (*,*) 'Treatment of periodic-data done'
290 
291  END SUBROUTINE prep_periodic_h_p_phi_bc
292 
293  SUBROUTINE list_periodic(np, jjs, sides, rr, side1, side2, e, list_out, perlist_out)
294  !============================================================================
295  IMPLICIT NONE
296  INTEGER, INTENT(IN) :: np
297  INTEGER, DIMENSION(:,:), INTENT(IN) :: jjs
298  INTEGER, DIMENSION(:), INTENT(IN) :: sides
299  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
300  INTEGER, INTENT(IN) :: side1, side2
301  REAL(KIND=8), DIMENSION(:), INTENT(IN) :: e
302  INTEGER, DIMENSION(:), POINTER :: list_out, perlist_out
303  INTEGER, DIMENSION(:), ALLOCATABLE :: list, perlist
304  LOGICAL, DIMENSION(np) :: virgin
305  REAL(KIND=8), DIMENSION(SIZE(rr,1)) :: ri
306  INTEGER :: ms, ns, i, j, long, inter
307  REAL(KIND=8) :: r, epsilon = 1.d-9
308  LOGICAL :: verif
309 
310  IF (ALLOCATED(list)) DEALLOCATE(list)
311  IF (ALLOCATED(perlist)) DEALLOCATE(perlist)
312 
313  ALLOCATE (list(np), perlist(np))
314  virgin = .true.
315 
316  i = 0; j=0
317  DO ms = 1, SIZE(sides)
318 
319  IF (sides(ms) .EQ. side1) THEN
320  DO ns = 1, SIZE(jjs,1)
321  IF (virgin(jjs(ns,ms))) THEN
322  i = i + 1
323  list(i) = jjs(ns,ms)
324  virgin(jjs(ns,ms)) = .false.
325  END IF
326  END DO
327  ELSE IF (sides(ms) .EQ. side2) THEN
328  DO ns = 1, SIZE(jjs,1)
329  IF (virgin(jjs(ns,ms))) THEN
330  j = j + 1
331  perlist(j) = jjs(ns,ms)
332  virgin(jjs(ns,ms)) = .false.
333  END IF
334  END DO
335 
336  END IF
337 
338  END DO
339 
340  IF (i .NE. j) THEN
341  WRITE(*,*) ' FEM_PERIODIC: side1 and side2 have', &
342  ' different numbers of points'
343  stop
344  END IF
345  long = i
346 
347  DO i = 1, long
348  ri = rr(:,list(i))+e(:)
349  verif = .false.
350  !if (i==2) stop
351  DO j = i, long
352  r = sum(abs(ri - rr(:,perlist(j))))
353  !if (i==1) write(*,*) ' r',r,'j', j
354  IF (r .LE. epsilon ) THEN
355  inter = perlist(i)
356  perlist(i) = perlist(j)
357  perlist(j) = inter
358  verif = .true.
359  EXIT
360  END IF
361  END DO
362  IF (.NOT.verif) THEN
363  WRITE(*,*) ' BUG dans data_periodic ou le maillage:', &
364  ' side1 + e /= side2'
365  WRITE(*,*) ' i = ', i
366  ! STOP
367  END IF
368  END DO
369 
370  ALLOCATE (list_out(long))
371  list_out(1:long) = list(1:long)
372  ALLOCATE (perlist_out(long))
373  perlist_out(1:long) = perlist(1:long)
374 
375  END SUBROUTINE list_periodic
376 
377  SUBROUTINE periodic_matrix_petsc(n_bord, list, perlist, matrix, LA)
378  USE dyn_line
379  USE def_type_mesh
380  USE my_util
381  IMPLICIT NONE
382  INTEGER , INTENT(IN) :: n_bord
383  TYPE(dyn_int_line), DIMENSION(:), INTENT(IN) :: list, perlist
384  TYPE(petsc_csr_la) , INTENT(IN) :: la
385  INTEGER, PARAMETER :: nmaxcols = 300
386  INTEGER :: ncols
387  INTEGER, DIMENSION(nmaxcols) :: cols
388  REAL(KIND=8), DIMENSION(nmaxcols) :: vals
389  INTEGER, DIMENSION(:), ALLOCATABLE :: n_cols_i
390  INTEGER, DIMENSION(1) :: idxn
391  INTEGER, DIMENSION(:,:), ALLOCATABLE :: jdxn
392  REAL(KIND=8), DIMENSION(:,:), ALLOCATABLE :: vals_pi
393  INTEGER :: n, l, i, pi, n_d, k
394 #include "petsc/finclude/petsc.h"
395  mat :: matrix
396  petscerrorcode :: ierr
397 
398  WRITE(*,*) 'Entering periodic_matrix_petsc'
399 
400  CALL matsetoption(matrix, mat_row_oriented, petsc_false, ierr)
401  CALL matsetoption(matrix, mat_keep_nonzero_pattern, petsc_true, ierr)
402 
403  DO k = 1, SIZE(la%loc_to_glob,1)
404  DO n = 1, n_bord
405  n_d = SIZE(list(n)%DIL)
406  IF (n_d /=0) THEN
407  ALLOCATE(jdxn(n_d,nmaxcols), vals_pi(n_d,nmaxcols), n_cols_i(n_d))
408  jdxn = 0
409  vals_pi = 0.d0
410  n_cols_i = 0
411 
412  DO l = 1, SIZE(list(n)%DIL)
413  idxn(1) = la%loc_to_glob(k,list(n)%DIL(l))
414  CALL matgetrow(matrix, idxn(1)-1,ncols,cols,vals,ierr)
415  n_cols_i(l) = ncols
416  jdxn(l,1:ncols) = cols(1:ncols)
417  vals_pi(l,1:ncols) = vals(1:ncols)
418  CALL matrestorerow(matrix, idxn(1)-1, ncols,cols, vals, ierr)
419  END DO
420 
421  DO l= 1, n_d
422  idxn(1) = la%loc_to_glob(k,perlist(n)%DIL(l))-1
423  CALL matsetvalues(matrix, 1,idxn, n_cols_i(l), jdxn(l,1:n_cols_i(l)), &
424  vals_pi(l:l,1:n_cols_i(l)), add_values, ierr)
425  END DO
426  DEALLOCATE(jdxn, vals_pi, n_cols_i)
427 
428  END IF
429  END DO
430  CALL matassemblybegin(matrix,mat_final_assembly,ierr)
431  CALL matassemblyend(matrix,mat_final_assembly,ierr)
432 
433  DO n = 1, n_bord
434  n_d = SIZE(list(n)%DIL)
435  CALL matzerorows(matrix, n_d, la%loc_to_glob(k,list(n)%DIL(:))-1, 1.d0, &
436  petsc_null_object, petsc_null_object, ierr)
437  END DO
438 !!$ CALL MatAssemblyBegin(matrix,MAT_FINAL_ASSEMBLY,ierr)
439 !!$ CALL MatAssemblyEnd(matrix,MAT_FINAL_ASSEMBLY,ierr)
440 
441  DO n = 1, n_bord
442  DO l = 1, SIZE(list(n)%DIL)
443  i = la%loc_to_glob(k,list(n)%DIL(l))
444  pi = la%loc_to_glob(k,perlist(n)%DIL(l))
445  CALL matsetvalue(matrix, i-1, pi-1, -1.d0, insert_values, ierr)
446  END DO
447  END DO
448  CALL matassemblybegin(matrix,mat_final_assembly,ierr)
449  CALL matassemblyend(matrix,mat_final_assembly,ierr)
450 
451  END DO
452 
453  END SUBROUTINE periodic_matrix_petsc
454 
455  SUBROUTINE periodic_rhs_petsc(n_bord, list, perlist, v_rhs, LA)
456  USE dyn_line
457  USE def_type_mesh
458  IMPLICIT NONE
459  INTEGER , INTENT(IN) :: n_bord
460  TYPE(dyn_int_line), DIMENSION(:), INTENT(IN) :: list, perlist
461  TYPE(petsc_csr_la) , INTENT(IN) :: la
462  INTEGER, DIMENSION(:), ALLOCATABLE :: idxn, jdxn
463  REAL(KIND=8), DIMENSION(:), ALLOCATABLE :: vals, bs
464  INTEGER :: n, k, n_d
465 #include "petsc/finclude/petsc.h"
466  vec :: v_rhs
467  petscerrorcode :: ierr
468 
469 
470  DO k = 1, SIZE(la%loc_to_glob,1)
471  DO n = 1, n_bord
472  n_d = SIZE(list(n)%DIL)
473  ALLOCATE(idxn(n_d), vals(n_d), jdxn(n_d), bs(n_d))
474  idxn = la%loc_to_glob(k,list(n)%DIL(:))-1
475  jdxn = la%loc_to_glob(k,perlist(n)%DIL(:))-1
476  CALL vecgetvalues(v_rhs, n_d, idxn, vals, ierr)
477  CALL vecassemblybegin(v_rhs,ierr)
478  CALL vecassemblyend(v_rhs,ierr)
479 
480  bs = 0.d0
481  CALL vecsetvalues(v_rhs, n_d, jdxn, vals, add_values, ierr)
482  CALL vecassemblybegin(v_rhs,ierr)
483  CALL vecassemblyend(v_rhs,ierr)
484  CALL vecsetvalues(v_rhs, n_d, idxn, bs, insert_values, ierr)
485  CALL vecassemblybegin(v_rhs,ierr)
486  CALL vecassemblyend(v_rhs,ierr)
487 
488  IF (ALLOCATED(idxn)) DEALLOCATE(idxn, jdxn, vals, bs)
489  END DO
490  END DO
491 
492  END SUBROUTINE periodic_rhs_petsc
493 
494 !!$ SUBROUTINE remove_extra_per(communicator, mesh, n_per_loc, list_loc)
495 !!$ USE def_type_mesh
496 !!$ IMPLICIT NONE
497 !!$ TYPE(mesh_type), INTENT(IN) :: mesh
498 !!$ INTEGER, INTENT(IN) :: n_per_loc
499 !!$ INTEGER, DIMENSION(n_per_loc), INTENT(INOUT) :: list_loc
500 !!$ INTEGER, DIMENSION(n_per_loc) :: node_loc
501 !!$ INTEGER, DIMENSION(:), ALLOCATABLE :: node_glob, list_n_per_loc, list_disp
502 !!$ INTEGER :: n_per_tot, rank, nb_proc, i
503 !!$#include "petsc/finclude/petsc.h"
504 !!$ MPI_Comm :: communicator
505 !!$ PetscErrorCode :: ierr
506 !!$
507 !!$ IF (n_per_loc /= 0) THEN
508 !!$ node_loc = mesh%loc_to_glob(list_loc)
509 !!$ END IF
510 !!$
511 !!$ CALL MPI_COMM_SIZE(communicator, nb_proc, ierr)
512 !!$ CALL MPI_COMM_RANK(communicator, rank, ierr)
513 !!$
514 !!$ ALLOCATE(list_n_per_loc(nb_proc),list_disp(nb_proc))
515 !!$
516 !!$ CALL MPI_ALLGATHER(n_per_loc, 1, MPI_INTEGER, list_n_per_loc, 1, MPI_INTEGER, communicator, ierr)
517 !!$ n_per_tot = SUM(list_n_per_loc)
518 !!$ ALLOCATE(node_glob(n_per_tot))
519 !!$
520 !!$ list_disp=0
521 !!$ DO i= 2,nb_proc
522 !!$ list_disp(i) = list_disp(i-1)+list_n_per_loc(i-1)
523 !!$ END DO
524 !!$
525 !!$ CALL MPI_ALLGATHERV(node_loc, n_per_loc, MPI_INTEGER, node_glob, list_n_per_loc, list_disp, MPI_INTEGER, &
526 !!$ communicator, ierr)
527 !!$
528 !!$ WRITE(200+rank,*) '---------------------------------'
529 !!$ WRITE(200+rank,*) node_glob
530 !!$ WRITE(200+rank,*) rank, node_loc
531 !!$ WRITE(200+rank,*) '---------------------------------'
532 !!$
533 !!$
534 !!$ DO i = 1, n_per_tot-1
535 !!$ IF (MINVAL(ABS(node_glob(i)-node_glob(i+1:))) == 0) THEN
536 !!$ node_glob(i) = -1
537 !!$ END IF
538 !!$ END DO
539 !!$
540 !!$ DO i = 1, n_per_loc
541 !!$ IF (node_glob(list_disp(rank+1)+i) == -1) THEN
542 !!$ list_loc(i) = -1
543 !!$ END IF
544 !!$ END DO
545 !!$
546 !!$ IF (ALLOCATED(node_glob)) DEALLOCATE(node_glob)
547 !!$ IF (ALLOCATED(list_n_per_loc)) DEALLOCATE(list_n_per_loc)
548 !!$ IF (ALLOCATED(list_disp)) DEALLOCATE(list_disp)
549 !!$
550 !!$ END SUBROUTINE remove_extra_per
551 
552 !!$ FUNCTION last_c_leng (len_str, string) RESULT (leng)
553 !!$ !===================================================
554 !!$ IMPLICIT NONE
555 !!$ INTEGER, INTENT(IN) :: len_str
556 !!$ CHARACTER (LEN=len_str), INTENT(IN) :: string
557 !!$ INTEGER :: leng
558 !!$ INTEGER :: i
559 !!$
560 !!$ leng = len_str
561 !!$
562 !!$ DO i=1,len_str
563 !!$ IF ( string(i:i) .EQ. ' ' ) THEN
564 !!$ leng = i-1; EXIT
565 !!$ ENDIF
566 !!$ ENDDO
567 !!$
568 !!$ END FUNCTION last_c_leng
569 
570 END MODULE periodic
subroutine list_periodic(np, jjs, sides, rr, side1, side2, e, list_out, perlist_out)
subroutine, public prep_periodic_bloc(my_periodic, mesh, periodic, nb_bloc)
subroutine, public periodic_rhs_petsc(n_bord, list, perlist, v_rhs, LA)
subroutine, public prep_periodic_h_p_phi_bc(my_periodic, H_mesh, pmag_mesh, phi_mesh, H_p_phi_per)
subroutine, public periodic_matrix_petsc(n_bord, list, perlist, matrix, LA)
subroutine, public prep_periodic(my_periodic, mesh, periodic)