OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
spmd_collect_nlocal.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!||====================================================================
24!|| spmd_collect_nlocal ../engine/source/mpi/output/spmd_collect_nlocal.F
25!||--- called by ------------------------------------------------------
26!|| resol ../engine/source/engine/resol.F
27!||--- calls -----------------------------------------------------
28!||--- uses -----------------------------------------------------
29!|| debug_mod ../engine/share/modules/debug_mod.F
30!|| inoutfile_mod ../common_source/modules/inoutfile_mod.F
31!|| nlocal_reg_mod ../common_source/modules/nlocal_reg_mod.F
32!|| spmd_comm_world_mod ../engine/source/mpi/spmd_comm_world.F90
33!||====================================================================
34 SUBROUTINE spmd_collect_nlocal(A ,SIZEA ,NUMNOD_LOCAL,
35 . POSI ,NLOC_DMG,SIZP0 ,
36 . NODGLOB,ITAB )
37
39 USE debug_mod
41C-----------------------------------------------
42C I m p l i c i t T y p e s
43C-----------------------------------------------
44 USE spmd_comm_world_mod, ONLY : spmd_comm_world
45#include "implicit_f.inc"
46C-----------------------------------------------------------------
47C M e s s a g e P a s s i n g
48C-----------------------------------------------
49#include "spmd.inc"
50C-----------------------------------------------
51C C o m m o n B l o c k s
52C-----------------------------------------------
53#include "com01_c.inc"
54#include "task_c.inc"
55#include "spmd_c.inc"
56#include "chara_c.inc"
57#include "units_c.inc"
58C-----------------------------------------------
59C D u m m y A r g u m e n t s
60C-----------------------------------------------
61 INTEGER ITAB(*),NUMNOD_LOCAL,NODGLOB(*)
62 INTEGER NPOS,SIZP0
63 INTEGER MSGOFF,MSGOFF0,MSGTYP,INFO,I,K,NG,N,
64 . empl,sdnodg(numnodm),filen
65 double precision
66 . aglob(2,numnodm),recglob(2,sizp0)
67C-----------------------------------------------
68C L o c a l V a r i a b l e s
69C-----------------------------------------------
70 INTEGER POSI(*)
71 INTEGER SIZEA
72 my_real
73 . a(sizea)
74 TYPE(NLOCAL_STR_), TARGET :: NLOC_DMG
75 CHARACTER FILNAM*100,CYCLENUM*7
76
77 INTEGER :: LEN_TMP_NAME
78 CHARACTER(len=2148) :: TMP_NAME
79
80#ifdef MPI
81 INTEGER STATUS(MPI_STATUS_SIZE),IERROR
82
83 DATA msgoff0/176/
84 DATA msgoff/177/
85C-----------------------------------------------
86C S o u r c e L i n e s
87C-----------------------------------------------
88 recglob(1:2,1:sizp0) = -1
89 IF (ispmd/=0) THEN
90 n=0
91 DO i = 1,numnod_local
92 npos = posi(i) ! Position of its first d.o.f in A
93
94 n=n+1
95 sdnodg(n) = nodglob(nloc_dmg%INDX(i))
96 aglob(1,n) = itab( nloc_dmg%INDX(i) )
97 aglob(2,n) = a(npos)
98 ENDDO
99
100
101 msgtyp=msgoff0
102 CALL mpi_send(sdnodg,n,mpi_integer,
103 . it_spmd(1),msgtyp,
104 . spmd_comm_world,ierror)
105 msgtyp=msgoff
106 CALL mpi_send(aglob,2*n,mpi_double_precision,
107 . it_spmd(1),msgtyp,
108 . spmd_comm_world,ierror)
109 ELSE
110
111 DO k=2,nspmd
112 msgtyp=msgoff0
113 CALL mpi_recv(sdnodg,numnodm,mpi_integer,
114 . it_spmd(k),msgtyp,
115 . spmd_comm_world,status,ierror)
116
117 CALL mpi_get_count(status,mpi_integer,n,ierror)
118
119 msgtyp=msgoff
120 CALL mpi_recv(aglob,2*n,mpi_double_precision,
121 . it_spmd(k),msgtyp,
122 . spmd_comm_world,status,ierror)
123
124
125 DO i=1,n
126 empl = sdnodg(i)
127 recglob(1,empl) = aglob(1,i)
128 recglob(2,empl) = aglob(2,i)
129 ENDDO
130 END DO
131 ENDIF
132#endif
133C
134
135 IF(ispmd==0) THEN
136 WRITE(cyclenum,'(I7.7)')ncycle
137 filnam=rootnam(1:rootlen)//'_NLOCAL_'//chrun//'_'//cyclenum//'.adb'
138
139 len_tmp_name = outfile_name_len + len_trim(filnam)
140 tmp_name=outfile_name(1:outfile_name_len)//filnam(1:len_trim(filnam))
141
142 OPEN(unit=idbg5,file=tmp_name(1:len_tmp_name),access='SEQUENTIAL',
143 . form='FORMATTED',status='UNKNOWN')
144
145 filen = rootlen+17
146
147 DO i=1,numnod_local
148 npos = posi(i) ! Position of its first d.o.f in A
149 n = nodglob(nloc_dmg%INDX(i))
150 recglob(1,n) = itab( nloc_dmg%INDX(i) )
151 recglob(2,n) = a(npos)
152 ENDDO
153
154
155 DO i = 1, numnodg
156 IF(nint(recglob(1,i))/=-1) THEN
157 WRITE(idbg5,'(A,I10,I10,Z20)' ) '>',ncycle,nint(recglob(1,i)),recglob(2,i)
158 ENDIF
159 END DO
160 WRITE (iout,1300) filnam(1:filen)
161 WRITE (istdo,1300) filnam(1:filen)
162 CLOSE(unit=idbg5)
163
164 END IF
165C
166 1300 FORMAT (4x,' DEBUG ANALYSIS NLOCAL FILE :',1x,a,' WRITTEN')
167 RETURN
168 END SUBROUTINE spmd_collect_nlocal
#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_collect_nlocal(a, sizea, numnod_local, posi, nloc_dmg, sizp0, nodglob, itab)