OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
debug_chksm_mod Module Reference

Functions/Subroutines

subroutine spmd_flush_accel (ncycle, ispmd, nspmd, numnod, numnodg, numnodm, a, itab, weight, nodglob)

Function/Subroutine Documentation

◆ spmd_flush_accel()

subroutine debug_chksm_mod::spmd_flush_accel ( integer, intent(in) ncycle,
integer, intent(in) ispmd,
integer, intent(in) nspmd,
integer, intent(in) numnod,
integer, intent(in) numnodg,
integer, intent(in) numnodm,
dimension(3,numnod), intent(in) a,
integer, dimension(numnod), intent(in) itab,
integer, dimension(numnod), intent(in) weight,
integer, dimension(numnod), intent(in) nodglob )
Parameters
[in]ncyclecurrent cycle
[in]ispmdMPI rank
[in]nspmdnumber of MPI processes
[in]numnodnumber of nodes on the current MPI process
[in]numnodgtotal number of nodes
[in]numnodmmaximum number of nodes over the MPI processes
[in]itabuser node id
[in]weightweight of the node (1 on the processor that assembles the forces, 0 elsewhere)
[in]nodglobglobal node id

Definition at line 42 of file spmd_flush_accel.F.

45C-----------------------------------------------
46C M o d u l e s
47C-----------------------------------------------
49C-----------------------------------------------
50C I m p l i c i t T y p e s
51C-----------------------------------------------
52 USE spmd_comm_world_mod, ONLY : spmd_comm_world
53#include "implicit_f.inc"
54C-----------------------------------------------
55C C o m m o n B l o c k s
56C-----------------------------------------------
57#include "units_c.inc"
58C-----------------------------------------------
59C M e s s a g e P a s s i n g
60C-----------------------------------------------
61#include "spmd.inc"
62C-----------------------------------------------
63C D u m m y A r g u m e n t s
64C-----------------------------------------------
65 INTEGER, INTENT(IN) :: NCYCLE !< current cycle
66 INTEGER, INTENT(IN) :: ISPMD !< MPI rank
67 INTEGER, INTENT(IN) :: NSPMD !< number of MPI processes
68 INTEGER, INTENT(IN) :: NUMNOD !< number of nodes on the current MPI process
69 INTEGER, INTENT(IN) :: NUMNODG !< total number of nodes
70 INTEGER, INTENT(IN) :: NUMNODM !< maximum number of nodes over the MPI processes
71 INTEGER, INTENT(IN) :: ITAB(NUMNOD) !< user node id
72 INTEGER, INTENT(IN) :: WEIGHT(NUMNOD) !< weight of the node (1 on the processor that assembles the forces, 0 elsewhere)
73 INTEGER, INTENT(IN) :: NODGLOB(NUMNOD) !< global node id
74 my_real, INTENT(IN) :: a(3,numnod) !< nodal acceleration
75C-----------------------------------------------
76C L o c a l V a r i a b l e s
77C-----------------------------------------------
78#ifdef MPI
79 INTEGER STATUS(MPI_STATUS_SIZE),IERROR
80#endif
81
82 INTEGER :: MSGOFF,MSGOFF0,MSGTYP !< MPI messages
83 INTEGER :: I,K,N,NODE_GLOBAL_ID !< loop indices
84 INTEGER, DIMENSION(:), ALLOCATABLE :: NODES_TO_SEND
85 DOUBLE PRECISION, DIMENSION(:,:), ALLOCATABLE :: AGLOB,NODES_TO_RECV
86 INTEGER :: CHECKSUM !< checksum of the nodal acceleration
87 DATA msgoff0/176/
88 DATA msgoff/177/
89C-----------------------------------------------
90C S o u r c e L i n e s
91C-----------------------------------------------
92
93 IF(ispmd == 0) THEN
94 ALLOCATE(nodes_to_recv(4,numnodg))
95 ENDIF
96
97 ALLOCATE(nodes_to_send(max(numnod,numnodm)))
98 ALLOCATE(aglob(4,max(numnodm,numnod)))
99
100 IF (ispmd/=0) THEN
101#ifdef MPI
102 n = 0
103 DO i = 1, numnod
104 IF (weight(i)==1) THEN
105 n = n+1
106 nodes_to_send(n) = nodglob(i)
107 aglob(1,n) = itab(i)
108 aglob(2,n) = a(1,i)
109 aglob(3,n) = a(2,i)
110 aglob(4,n) = a(3,i)
111 END IF
112 END DO
113 msgtyp=msgoff0
114 CALL mpi_send(nodes_to_send,n,mpi_integer,
115 . 0,msgtyp,
116 . spmd_comm_world,ierror)
117 msgtyp=msgoff
118 CALL mpi_send(aglob,4*n,mpi_double_precision,
119 . 0,msgtyp,
120 . spmd_comm_world,ierror)
121#endif
122 ELSE ! ISPMD == 0
123
124 DO i=1,numnod
125 IF (weight(i)==1) THEN
126 node_global_id = nodglob(i)
127 nodes_to_recv(1,node_global_id) = itab(i)
128 nodes_to_recv(2,node_global_id) = a(1,i)
129 nodes_to_recv(3,node_global_id) = a(2,i)
130 nodes_to_recv(4,node_global_id) = a(3,i)
131 ENDIF
132 ENDDO
133#ifdef MPI
134 DO k=2,nspmd
135 msgtyp=msgoff0
136 CALL mpi_recv(nodes_to_send,numnodm,mpi_integer,
137 . k-1,msgtyp,
138 . spmd_comm_world,status,ierror)
139 CALL mpi_get_count(status,mpi_integer,n,ierror)
140 msgtyp=msgoff
141 CALL mpi_recv(aglob,4*n,mpi_double_precision,
142 . k-1,msgtyp,
143 . spmd_comm_world,status,ierror)
144 DO i=1,n
145 node_global_id = nodes_to_send(i)
146 nodes_to_recv(1,node_global_id) = aglob(1,i)
147 nodes_to_recv(2,node_global_id) = aglob(2,i)
148 nodes_to_recv(3,node_global_id) = aglob(3,i)
149 nodes_to_recv(4,node_global_id) = aglob(4,i)
150 ENDDO
151
152 END DO
153#endif
154 checksum = double_array_checksum(nodes_to_recv,numnodg,4)
155 WRITE(iout,*) ncycle, "CHECKSUM:",checksum
156
157 ENDIF
158 IF(ALLOCATED(nodes_to_send)) DEALLOCATE(nodes_to_send)
159 IF(ALLOCATED(aglob)) DEALLOCATE(aglob)
160 IF(ALLOCATED(nodes_to_recv)) DEALLOCATE(nodes_to_recv)
161 RETURN
#define my_real
Definition cppsort.cpp:32
#define max(a, b)
Definition macros.h:21
subroutine mpi_recv(buf, cnt, datatype, source, tag, comm, status, ierr)
Definition mpi.f:461
subroutine mpi_get_count(status, datatype, cnt, ierr)
Definition mpi.f:296
subroutine mpi_send(buf, cnt, datatype, dest, tag, comm, ierr)
Definition mpi.f:480
integer function double_array_checksum(a, siz1, siz2)