OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
spmd_collectt.F
Go to the documentation of this file.
1Copyright> OpenRadioss
2Copyright> Copyright (C) 1986-2025 Altair Engineering Inc.
3Copyright>
4Copyright> This program is free software: you can redistribute it and/or modify
5Copyright> it under the terms of the GNU Affero General Public License as published by
6Copyright> the Free Software Foundation, either version 3 of the License, or
7Copyright> (at your option) any later version.
8Copyright>
9Copyright> This program is distributed in the hope that it will be useful,
10Copyright> but WITHOUT ANY WARRANTY; without even the implied warranty of
11Copyright> MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12Copyright> GNU Affero General Public License for more details.
13Copyright>
14Copyright> You should have received a copy of the GNU Affero General Public License
15Copyright> along with this program. If not, see <https://www.gnu.org/licenses/>.
16Copyright>
17Copyright>
18Copyright> Commercial Alternative: Altair Radioss Software
19Copyright>
20Copyright> As an alternative to this open-source version, Altair also offers Altair Radioss
21Copyright> software under a commercial license. Contact Altair to discuss further if the
22Copyright> commercial version may interest you: https://www.altair.com/radioss/.
23
24C
25!||====================================================================
26!|| spmd_collectt ../engine/source/mpi/output/spmd_collectt.F
27!||--- called by ------------------------------------------------------
28!|| resol ../engine/source/engine/resol.F
29!||--- calls -----------------------------------------------------
30!||--- uses -----------------------------------------------------
31!|| inoutfile_mod ../common_source/modules/inoutfile_mod.F
32!|| spmd_comm_world_mod ../engine/source/mpi/spmd_comm_world.F90
33!||====================================================================
34 SUBROUTINE spmd_collectt(TEMP, ITAB, WEIGHT,NODGLOB,SIZP0)
36C-----------------------------------------------
37C I m p l i c i t T y p e s
38C-----------------------------------------------
39 USE spmd_comm_world_mod, ONLY : spmd_comm_world
40#include "implicit_f.inc"
41C-----------------------------------------------------------------
42C M e s s a g e P a s s i n g
43C-----------------------------------------------
44#include "spmd.inc"
45C-----------------------------------------------
46C C o m m o n B l o c k s
47C-----------------------------------------------
48#include "com01_c.inc"
49#include "com04_c.inc"
50#include "task_c.inc"
51#include "spmd_c.inc"
52#include "chara_c.inc"
53#include "units_c.inc"
54C-----------------------------------------------
55C D u m m y A r g u m e n t s
56C-----------------------------------------------
57 INTEGER ITAB(*),WEIGHT(*),NODGLOB(*),SIZP0
59 . temp(*)
60C-----------------------------------------------
61C L o c a l V a r i a b l e s
62C-----------------------------------------------
63#ifdef MPI
64 INTEGER MSGOFF,MSGOFF0,MSGTYP,INFO,I,K,NG,N,
65 . EMPL,SDNODG(NUMNODM),FILEN
66 INTEGER STATUS(MPI_STATUS_SIZE),IERROR
67 double precision
68 . aglob(2,numnodm),recglob(2,sizp0)
69
70 CHARACTER FILNAM*100,CYCLENUM*7
71
72 INTEGER :: LEN_TMP_NAME
73 CHARACTER(len=2148) :: TMP_NAME
74
75 DATA msgoff0/184/
76 DATA msgoff/185/
77C-----------------------------------------------
78C S o u r c e L i n e s
79C-----------------------------------------------
80 WRITE(cyclenum,'(I7.7)')ncycle
81 filnam=rootnam(1:rootlen)//'_'//chrun//'_'//cyclenum//'.tdb'
82
83 len_tmp_name = outfile_name_len + len_trim(filnam)
84 tmp_name=outfile_name(1:outfile_name_len)//filnam(1:len_trim(filnam))
85
86 OPEN(unit=idbg8,file=tmp_name(1:len_tmp_name),access='SEQUENTIAL',
87 . form='FORMATTED',status='UNKNOWN')
88
89 filen = rootlen+17
90 IF (ispmd/=0) THEN
91 n = 0
92 DO i = 1, numnod
93 IF (weight(i)==1) THEN
94 n = n+1
95 sdnodg(n) = nodglob(i)
96 aglob(1,n) = itab(i)
97 aglob(2,n) = temp(i)
98 END IF
99 END DO
100 msgtyp=msgoff0
101 CALL mpi_send(sdnodg,n,mpi_integer,
102 . it_spmd(1),msgtyp,
103 . spmd_comm_world,ierror)
104 msgtyp=msgoff
105 CALL mpi_send(aglob,2*n,mpi_double_precision,
106 . it_spmd(1),msgtyp,
107 . spmd_comm_world,ierror)
108 ELSE
109 DO i=1,numnod
110 IF (weight(i)==1) THEN
111 empl = nodglob(i)
112 recglob(1,empl) = itab(i)
113 recglob(2,empl) = temp(i)
114 ENDIF
115 ENDDO
116
117 DO k=2,nspmd
118 msgtyp=msgoff0
119 CALL mpi_recv(sdnodg,numnodm,mpi_integer,
120 . it_spmd(k),msgtyp,
121 . spmd_comm_world,status,ierror)
122
123 CALL mpi_get_count(status,mpi_integer,n,ierror)
124
125 msgtyp=msgoff
126 CALL mpi_recv(aglob,2*n,mpi_double_precision,
127 . it_spmd(k),msgtyp,
128 . spmd_comm_world,status,ierror)
129
130 DO i=1,n
131 empl = sdnodg(i)
132 recglob(1,empl) = aglob(1,i)
133 recglob(2,empl) = aglob(2,i)
134 ENDDO
135
136 END DO
137C
138 DO i = 1, numnodg
139 WRITE(idbg8,'(A,I10,I10,Z20)' )
140 . '>',ncycle,nint(recglob(1,i)),
141 . recglob(2,i)
142 END DO
143 WRITE (iout,1300) filnam(1:filen)
144 WRITE (istdo,1300) filnam(1:filen)
145 CLOSE(unit=idbg8)
146
147 END IF
148C
149 1300 FORMAT (4x,' DEBUG ANALYSIS FILE:',1x,a,' WRITTEN')
150#endif
151 RETURN
152 END
#define my_real
Definition cppsort.cpp:32
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
character(len=outfile_char_len) outfile_name
integer outfile_name_len
subroutine spmd_collectt(temp, itab, weight, nodglob, sizp0)