39#include "Flashx_mpi_implicitNone.fh"
42 integer,
intent(in) :: myPE
48 integer :: istat, i, strtype
49 integer :: n_real, n_int, n_str, n_log
50 real,
allocatable :: real_vals(:)
51 integer,
allocatable :: int_vals(:)
52 logical,
allocatable :: log_vals(:)
53 character(
len=MAX_STRING_LENGTH),
allocatable :: str_vals(:),
&
62 call MPI_Type_Contiguous (MAX_STRING_LENGTH,
MPI_CHARACTER, strtype, istat)
63 call MPI_Type_Commit (strtype, istat)
74 n_real
= context
%n_real
75 call MPI_Bcast (n_real,
1,
MPI_INTEGER,
MASTER_PE,
&
76 & MPI_COMM_WORLD, istat)
80 allocate (real_names(n_real), real_vals(n_real),
stat=istat)
82 write (
*,
*)
'nameValueLL_bcast : allocate failed in real bcast'
83 call Driver_abort(
"Error: nameValueLL_bcast : allocate failed in real bcast")
88 node_real
=> context
%real_list
89 do while (
associated(node_real))
90 real_names(i)
= node_real
%name
91 real_vals(i)
= node_real
%value
93 node_real
=> node_real
%next
97 call MPI_Bcast (real_names, n_real,
&
98 & strtype,
MASTER_PE,
MPI_COMM_WORLD, istat)
99 call MPI_Bcast (real_vals, n_real,
&
101 & MPI_COMM_WORLD, istat)
109 deallocate (real_names, real_vals)
116 n_int
= context
%n_int
117 call MPI_Bcast (n_int,
1,
MPI_INTEGER,
MASTER_PE,
&
118 & MPI_COMM_WORLD, istat)
122 allocate (int_names(n_int), int_vals(n_int),
stat=istat)
124 write (
*,
*)
'nameValueLL_bcast: allocate failed'
125 call Driver_abort(
"Error: nameValueLL_bcast : allocate failed")
130 node_int
=> context
%int_list
131 do while (
associated(node_int))
132 int_names(i)
= node_int
%name
133 int_vals(i)
= node_int
%value
135 node_int
=> node_int
%next
139 call MPI_Bcast (int_names, n_int,
&
140 & strtype,
MASTER_PE,
MPI_COMM_WORLD, istat)
141 call MPI_Bcast (int_vals, n_int,
&
143 & MPI_COMM_WORLD, istat)
151 deallocate (int_names, int_vals)
157 n_str
= context
%n_str
158 call MPI_Bcast (n_str,
1,
MPI_INTEGER,
MASTER_PE,
&
159 & MPI_COMM_WORLD, istat)
163 allocate (str_names(n_str), str_vals(n_str),
stat=istat)
165 write (
*,
*)
'nameValueLL_bcast : allocate failed'
166 call Driver_abort(
"Error: nameValueLL_bcast : allocate failed");
171 node_str
=> context
%str_list
172 do while (
associated(node_str))
173 str_names(i)
= node_str
%name
174 str_vals(i)
= node_str
%value
176 node_str
=> node_str
%next
180 call MPI_Bcast (str_names, n_str,
&
181 & strtype,
MASTER_PE,
MPI_COMM_WORLD, istat)
182 call MPI_Bcast (str_vals, n_str,
&
183 & strtype,
MASTER_PE,
MPI_COMM_WORLD, istat)
191 deallocate (str_names, str_vals)
197 n_log
= context
%n_log
198 call MPI_Bcast (n_log,
1,
MPI_INTEGER,
MASTER_PE,
&
199 & MPI_COMM_WORLD, istat)
203 allocate (log_names(n_log), log_vals(n_log),
stat=istat)
205 write (
*,
*)
'nameValueLL_bcast : allocate failed'
206 call Driver_abort(
"nameValueLL_bcast : allocate failed")
211 node_log
=> context
%log_list
212 do while (
associated(node_log))
213 log_names(i)
= node_log
%name
214 log_vals(i)
= node_log
%value
216 node_log
=> node_log
%next
220 call MPI_Bcast (log_names, n_log,
&
221 & strtype,
MASTER_PE,
MPI_COMM_WORLD, istat)
222 call MPI_Bcast (log_vals, n_log,
&
224 & MPI_COMM_WORLD, istat)
232 deallocate (log_names, log_vals)
238 call MPI_Type_Free (strtype, istat)
subroutine nameValueLL_bcast(context, myPE)
integer, parameter TYPE_VAR