progam.f90
1 ! $Id: progam.f90,v 1.2 2006-07-17 00:52:06 emiliano Exp $
2 ! PROGAM-P - Programa de Generacin Automática de Mallas en Paralelo 3 4 PROGRAM progam_p 5 !DEC$ REAL:8 6 USE region_mod 7 IMPLICIT NONE 8 9 TYPE(region_t) :: region_local 10 INTEGER :: ierr 11 12 CALL MPI_INIT(ierr) 13
14 CALL region_init_mod() ! Inicialización 15 CALL region_input_p(region_local) ! Entrada de datos
16 CALL region_grid(region_local) ! Creación de los nodos interiores 17 CALL region_nodenums_p(region_local) ! Numeración de los nodos
18 CALL region_elements(region_local) ! Creación de los elementos 19 CALL region_output_p(region_local) ! Salida para GiD
20
21 CALL MPI_FINALIZE(ierr) 22
23 END PROGRAM progam_p
region.f90
1 ! $Id: region.f90,v 1.2 2006-07-17 00:52:06 emiliano Exp $ 2 ! PROGAM-P(Programa de Generacin Automática de Mallas en Paralelo) 3 4 MODULE region_mod 5 !DEC$ REAL:8 6 IMPLICIT NONE 7 INCLUDE 'mpif.h' 8 9 TYPE region_t
10 INTEGER :: num ! número de región
11 INTEGER :: n, m ! cantidad de filas y columnas 12 INTEGER :: elem_type ! 3=triángulo, 4=cuadrado 13 INTEGER :: vecinas(4) ! datos de conectividad
14 REAL :: def_nodes(8,2) ! los 8 nodos que definen la región 15 REAL, ALLOCATABLE :: node_coords(:,:,:) ! coordenadas de los nodos
16 INTEGER, ALLOCATABLE :: node_nums(:,:) ! números de los nodos
17 INTEGER, ALLOCATABLE :: elements(:,:) ! cada elemento es un vector de 3 o 4 enteros 18 END TYPE region_t
19
20 INTEGER, PRIVATE :: my_rank, p, status(MPI_STATUS_SIZE), ierr 21
22 CONTAINS
23 ! Inicialización
24 SUBROUTINE region_init_mod()
25 CALL MPI_COMM_RANK(MPI_COMM_WORLD, my_rank, ierr) 26 CALL MPI_COMM_SIZE(MPI_COMM_WORLD, p, ierr) 27 END SUBROUTINE region_init_mod
28
29 ! Procesamiento de la entrada de datos 30 SUBROUTINE region_input_p(this) 31 TYPE(region_t), INTENT(OUT) :: this 32 INTEGER :: i
36 INTEGER :: vecinas(4) 37 INTEGER :: node_nums(8) 38
39 IF (my_rank == 0) THEN
40 READ(*,*) n_regions, n_input_nodes 41 IF (n_regions /= p) THEN
42 WRITE(*,*) 'El número de procesos es distinto del número de regiones.' 43 CALL MPI_ABORT(MPI_COMM_WORLD, 1, ierr)
44 ENDIF 45 ENDIF 46
47 CALL MPI_BCAST(n_regions, 1, MPI_INTEGER, 0, MPI_COMM_WORLD, ierr) 48 CALL MPI_BCAST(n_input_nodes, 1, MPI_INTEGER, 0, MPI_COMM_WORLD, ierr) 49 ALLOCATE(input_nodes(n_input_nodes,2))
50
51 ! Lectura de los nodos y sus coordenadas 52 IF (my_rank == 0) THEN
53 DO i=1, n_input_nodes
54 READ(*,*) node_num, input_nodes(i,:) 55 END DO
56 END IF 57
58 ! Los nodos se envían a todos los procesos, sin importar cuáles necesiten. 59 CALL MPI_BCAST(input_nodes, SIZE(input_nodes), MPI_REAL8, 0, MPI_COMM_WORLD, ierr) 60
61 IF (my_rank == 0) THEN
62 ! El proceso 0 se queda con la primera región.
63 READ(*,*) this%n, this%m, this%elem_type, this%vecinas, node_nums 64 this%def_nodes = input_nodes(node_nums,:)
65
66 ! Lee los datos de las regiones restantes y los envía al proceso 67 ! correspondiente.
68 DO i=2, n_regions
69 READ(*,*) n, m, elem_type, vecinas, node_nums
70 CALL MPI_SEND(n, 1, MPI_INTEGER, i-1, 0, MPI_COMM_WORLD, ierr) 71 CALL MPI_SEND(m, 1, MPI_INTEGER, i-1, 0, MPI_COMM_WORLD, ierr) 72 CALL MPI_SEND(elem_type, 1, MPI_INTEGER, i-1, 0, MPI_COMM_WORLD, ierr) 73 CALL MPI_SEND(vecinas, 4, MPI_INTEGER, i-1, 0, MPI_COMM_WORLD, ierr) 74 CALL MPI_SEND(node_nums, 8, MPI_INTEGER, i-1, 0, MPI_COMM_WORLD, ierr) 75 END DO
76
77 ELSE
78 ! Recibe los datos enviados por el proceso 0.
79 CALL MPI_RECV(this%n, 1, MPI_INTEGER, 0, 0, MPI_COMM_WORLD, status, ierr) 80 CALL MPI_RECV(this%m, 1, MPI_INTEGER, 0, 0, MPI_COMM_WORLD, status, ierr)
81 CALL MPI_RECV(this%elem_type, 1, MPI_INTEGER, 0, 0, MPI_COMM_WORLD, status, ierr) 82 CALL MPI_RECV(this%vecinas, 4, MPI_INTEGER, 0, 0, MPI_COMM_WORLD, status, ierr) 83 CALL MPI_RECV(node_nums, 8, MPI_INTEGER, 0, 0, MPI_COMM_WORLD, status, ierr) 84 this%def_nodes = input_nodes(node_nums,:)
85 END IF 86
87 this%num = my_rank + 1 88
89 END SUBROUTINE region_input_p 90
91 ! Divide la región en (n-1)*(m-1) elementos cuadrangulares. 92 SUBROUTINE region_grid(this)
93 TYPE(region_t), INTENT(INOUT) :: this 94 INTEGER :: i, j
95 REAL :: pun(8) ! funciones de forma 96 REAL :: eta, si ! coordenadas normalizadas 97 REAL :: deta, dsi ! incrementos de eta y si 98 99 ALLOCATE(this%node_coords(this%n,this%m,2)) 100 ALLOCATE(this%node_nums(this%n,this%m)) 101 102 deta = 2./(this%n-1) 103 dsi = 2./(this%m-1) 104 DO i=1, this%n 105 eta = 1 - (i-1)*deta 106 DO j=1, this%m 107 si = -1 + (j-1)*dsi
108 pun(1) = -0.25 * (1-si) * (1-eta) * (si+eta+1) 109 pun(2) = 0.5 * (1-si**2) * (1-eta)
110 pun(3) = 0.25 * (1+si) * (1-eta) * (si-eta-1) 111 pun(4) = 0.5 * (1+si) * (1-eta**2)
112 pun(5) = 0.25 * (1+si) * (1+eta) * (si+eta-1) 113 pun(6) = 0.5 * (1-si**2) * (1+eta)
114 pun(7) = 0.25 * (1-si) * (1+eta) * (eta-si-1) 115 pun(8) = 0.5 * (1-si) * (1 - eta**2)
116 this%node_coords(i,j,1) = SUM(this%def_nodes(:,1) * pun) 117 this%node_coords(i,j,2) = SUM(this%def_nodes(:,2) * pun) 118 END DO
119 END DO
120 END SUBROUTINE region_grid 121
122 ! Numera en forma global los nodos de la región, teniendo en cuenta los que 123 ! ya han sido numerados por otras regiones.
124 SUBROUTINE region_nodenums_p(this)
125 TYPE(region_t), INTENT(INOUT), TARGET :: this 126 INTEGER :: i, j, istart, jstart, iend, jend
127 INTEGER :: nstart ! La numeración comienza a partir de este número. 128 INTEGER :: nend ! Este valor se envía a la siguiente región. 129 INTEGER :: ranks(4) ! Ranks de las regiones vecinas.
130 INTEGER, POINTER :: bound_nums(:)
131 REAL :: chk_node(2), recv_node(2) 132 133 nstart = 1 134 istart = 1 135 jstart = 1 136 iend = this%n 137 jend = this%m 138
139 ! Verifica si los nodos de las fronteras son numerados por otra región. 140 IF (this%vecinas(1)<this%num .AND. this%vecinas(1)/=0) iend=iend-1 141 IF (this%vecinas(2)<this%num .AND. this%vecinas(2)/=0) jend=jend-1 142 IF (this%vecinas(3)<this%num .AND. this%vecinas(3)/=0) istart=istart+1 143 IF (this%vecinas(4)<this%num .AND. this%vecinas(4)/=0) jstart=jstart+1 144
145 ! Obtiene el valor nstart
146 IF (my_rank>0) CALL MPI_RECV(nstart, 1, MPI_INTEGER, my_rank-1, 0, MPI_COMM_WORLD, status, ierr)
147 nend = nstart + (iend-istart+1) * (jend-jstart+1)
148 IF (my_rank<p-1) CALL MPI_SEND(nend, 1, MPI_INTEGER, my_rank+1, 0, MPI_COMM_WORLD, ierr) 149
150 ! Numera los nodos 151 DO i=istart, iend 152 DO j=jstart, jend 153 this%node_nums(i,j) = nstart 154 nstart = nstart+1 155 END DO 156 END DO 157
158 ! Recibe los nodos de las fronteras numeradas por una región vecina. 159 ! El chk_node es para verificar que el vector recibido esté en el orden 160 ! correcto. Los mpi_recv se hacen antes que los mpi_send para contemplar 161 ! el caso de que un nodo ubicado en una esquina, sea numerado por una 162 ! región vecina y deba ser enviado a otra región vecina.
163 ranks = this%vecinas - 1 164 DO i=1, 4
165 SELECT CASE(i) 166 CASE(1)
167 bound_nums => this%node_nums(this%n,:) ! abajo 168 chk_node = this%node_coords(this%n,1,:)
169 CASE(2)
170 bound_nums => this%node_nums(:,this%m) ! derecha 171 chk_node = this%node_coords(1,this%m,:)
172 CASE(3)
173 bound_nums => this%node_nums(1,:) ! arriba 174 chk_node = this%node_coords(1,1,:)
175 CASE(4)
179
180 IF (this%vecinas(i)<this%num .AND. this%vecinas(i)/=0) THEN
181 CALL MPI_RECV(bound_nums, SIZE(bound_nums), MPI_INTEGER, ranks(i), 0, MPI_COMM_WORLD, status, ierr)
182 CALL MPI_RECV(recv_node, 2, MPI_REAL8, ranks(i), 1, MPI_COMM_WORLD, status, ierr) 183 IF (ANY(chk_node/=recv_node)) CALL invertir(bound_nums)
184 END IF 185 END DO 186
187 ! Envía los números de las fronteras a las regiones vecinas que tengan mayor 188 ! número de región.
189 DO i=1, 4 190 SELECT CASE(i) 191 CASE(1)
192 bound_nums => this%node_nums(this%n,:) ! abajo 193 chk_node = this%node_coords(this%n,1,:)
194 CASE(2)
195 bound_nums => this%node_nums(:,this%m) ! derecha 196 chk_node = this%node_coords(1,this%m,:)
197 CASE(3)
198 bound_nums => this%node_nums(1,:) ! arriba 199 chk_node = this%node_coords(1,1,:)
200 CASE(4)
201 bound_nums => this%node_nums(:,1) ! izquierda 202 chk_node = this%node_coords(1,1,:)
203 END SELECT 204
205 IF (this%vecinas(i)>this%num) THEN
206 CALL MPI_SEND(bound_nums, SIZE(bound_nums), MPI_INTEGER, ranks(i), 0, MPI_COMM_WORLD, ierr)
207 CALL MPI_SEND(chk_node, 2, MPI_REAL8, ranks(i), 1, MPI_COMM_WORLD, ierr) 208 END IF
209 END DO
210 END SUBROUTINE region_nodenums_p 211
212 ! Invierte el orden de los elementos en el vector. 213 SUBROUTINE invertir(vector)
214 INTEGER, INTENT(INOUT) :: vector(:) 215 INTEGER :: aux(SIZE(vector)) 216 INTEGER :: i, n 217 218 n = SIZE(vector) 219 DO i=1, n 220 aux(n-i+1) = vector(i) 221 END DO 222 vector = aux
223 END SUBROUTINE invertir 224
225 ! Genera elementos triangulares (dividiendo los cuadrados según la diagonal más 226 ! corta del elemento (1,1)) o cuadrangulares. Cada elemento es un vector formado 227 ! por los números de los nodos que lo definen.
228 SUBROUTINE region_elements(this)
229 TYPE(region_t), INTENT(INOUT) :: this 230 INTEGER :: i, j, k 231 REAL :: diag1, diag2 232
233 IF (this%elem_type == 3) THEN ! elementos triangulares 234 ALLOCATE(this%elements((this%n-1)*(this%m-1)*2, 3)) 235 k = 1
236 diag1 = SQRT((this%node_coords(1,1,1)-this%node_coords(2,2,1))**2 & 237 + (this%node_coords(1,1,2)-this%node_coords(2,2,2))**2) 238 diag2 = SQRT((this%node_coords(2,1,1)-this%node_coords(1,2,1))**2 & 239 + (this%node_coords(2,1,2)-this%node_coords(1,2,2))**2) 240
241 IF (diag1 < diag2) THEN 242 DO i=1, this%n-1 243 DO j=1, this%m-1 244 this%elements(k,1) = this%node_nums(i,j) 245 this%elements(k,2) = this%node_nums(i+1,j) 246 this%elements(k,3) = this%node_nums(i+1,j+1) 247 k = k+1 248 this%elements(k,1) = this%node_nums(i,j)
249 this%elements(k,2) = this%node_nums(i+1,j+1) 250 this%elements(k,3) = this%node_nums(i,j+1) 251 k = k+1 252 END DO 253 END DO 254 ELSE 255 DO i=1, this%n-1 256 DO j=1, this%m-1 257 this%elements(k,1) = this%node_nums(i,j) 258 this%elements(k,2) = this%node_nums(i+1,j) 259 this%elements(k,3) = this%node_nums(i,j+1) 260 k = k+1 261 this%elements(k,1) = this%node_nums(i+1,j) 262 this%elements(k,2) = this%node_nums(i+1,j+1) 263 this%elements(k,3) = this%node_nums(i,j+1) 264 k = k+1 265 END DO 266 END DO 267 END IF 268
269 ELSE IF (this%elem_type == 4) THEN ! elementos cuadrangulares 270 ALLOCATE(this%elements((this%n-1)*(this%m-1), 4)) 271 k = 1 272 DO i=1, this%n-1 273 DO j=1, this%m-1 274 this%elements(k,1) = this%node_nums(i,j) 275 this%elements(k,2) = this%node_nums(i+1,j) 276 this%elements(k,3) = this%node_nums(i+1,j+1) 277 this%elements(k,4) = this%node_nums(i,j+1) 278 k = k+1 279 END DO 280 END DO 281 END IF
282 END SUBROUTINE region_elements 283
284 ! Salida para GiD
285 SUBROUTINE region_output_p(this) 286 TYPE(region_t), INTENT(IN) :: this
287 INTEGER :: i, j, elem_num, node_num 288
289 IF (my_rank>0) THEN
290 CALL MPI_RECV(elem_num, 1, MPI_INTEGER, my_rank-1, 0, MPI_COMM_WORLD, status, ierr) 291 CALL MPI_RECV(node_num, 1, MPI_INTEGER, my_rank-1, 0, MPI_COMM_WORLD, status, ierr) 292 OPEN(UNIT=1, FILE="output.msh", ACTION='WRITE', POSITION='APPEND', STATUS='OLD') 293 ELSE
294 elem_num = 1 295 node_num = 0
296 OPEN(UNIT=1, FILE="output.msh", ACTION='WRITE', STATUS='REPLACE') 297 END IF
298 299
300 ! Cabecera
301 IF (this%elem_type==3) THEN
302 WRITE(1,*) 'MESH dimension 2 ElemType Triangle Nnode 3' 303 ELSE IF (this%elem_type==4) THEN
304 WRITE(1,*) 'MESH dimension 2 ElemType Quadrilateral Nnode 4' 305 END IF
306
307 ! Coordenadas de los nodos 308 WRITE(1,*) 'coordinates' 309 WRITE(1,*) '# num x y' 310 DO i=1, this%n 311 DO j=1, this%m 312 IF (this%node_nums(i,j)>node_num) THEN 313 node_num = this%node_nums(i,j)
314 WRITE(1,*) node_num, this%node_coords(i,j,:) 315 END IF
316 END DO 317 END DO
321 WRITE(1,*) 'elements'
322 WRITE(1,*) '# num nodo1 nodo2 nodo3 nodo4' 323 DO i=1, SIZE(this%elements, DIM=1)
324 WRITE(1,*) elem_num, this%elements(i,:) 325 elem_num = elem_num+1
326 END DO
327 WRITE(1,*) 'end elements' 328 CLOSE(1)
329
330 IF (my_rank<p-1) THEN
331 CALL MPI_SEND(elem_num, 1, MPI_INTEGER, my_rank+1, 0, MPI_COMM_WORLD, ierr) 332 CALL MPI_SEND(node_num, 1, MPI_INTEGER, my_rank+1, 0, MPI_COMM_WORLD, ierr) 333 END IF
334 END SUBROUTINE region_output_p 335 END MODULE region_mod
test.f90
1 ! $Id: test.f90,v 1.2 2006-07-17 00:52:06 emiliano Exp $
2 ! PROGAM-P - Programa de Generacin Automática de Mallas en Paralelo 3 4 PROGRAM progam_p_test 5 !DEC$ REAL:8 6 USE region_mod 7 IMPLICIT NONE 8 9 TYPE(region_t) :: region_local 10 INTEGER :: my_rank, p, ierr 11 REAL :: ti, tf
12 INTEGER :: n_local, m_local ! Valores iniciales de n y m 13 INTEGER :: x ! Factor de multiplicación de n y m 14 INTEGER, ALLOCATABLE :: n(:), m(:) ! Para calcular el número total de nodos 15
16 CALL MPI_INIT(ierr)
17 CALL MPI_COMM_RANK(MPI_COMM_WORLD, my_rank, ierr) 18 CALL MPI_COMM_SIZE(MPI_COMM_WORLD, p, ierr) 19
20 CALL region_init_mod() ! Inicialización 21 CALL region_input_p(region_local) ! Entrada de datos 22
23 n_local = region_local%n 24 m_local = region_local%m
25 IF (my_rank==0) ALLOCATE(n(p), m(p))
26 CALL MPI_GATHER(n_local, 1, MPI_INTEGER, n, 1, MPI_INTEGER, 0, MPI_COMM_WORLD, ierr) 27 CALL MPI_GATHER(m_local, 1, MPI_INTEGER, m, 1, MPI_INTEGER, 0, MPI_COMM_WORLD, ierr) 28
29 DO x=1, 100
30 region_local%n = n_local*x 31 region_local%m = m_local*x 32
33 CALL MPI_BARRIER(MPI_COMM_WORLD, ierr) 34 ti = MPI_WTIME()
35 CALL region_grid(region_local) ! Creación de los nodos interiores 36 CALL region_nodenums_p(region_local) ! Numeración de los nodos
37 CALL region_elements(region_local) ! Creación de los elementos 38 CALL MPI_BARRIER(MPI_COMM_WORLD, ierr)
39 tf = MPI_WTIME() 40
41 IF (my_rank==0) WRITE(*,*) DOT_PRODUCT(n*x, m*x), tf-ti 42 DEALLOCATE(region_local%node_nums) 43 DEALLOCATE(region_local%node_coords) 44 DEALLOCATE(region_local%elements) 45 END DO 46 47 CALL MPI_FINALIZE(ierr) 48