Status: | Merged |
---|---|
Merged at revision: | 1823caa92a161268c23f9b1379010796504c2231 |
Proposed branch: | libfdf:ag |
Merge into: | libfdf:master |
Diff against target: |
1320 lines (+246/-669) 5 files modified
README.md (+8/-7) configure.ac (+1/-1) doc/fdf.Standard (+0/-17) src/fdf.F90 (+235/-626) src/utils.F90 (+2/-18) |
Related bugs: |
Reviewer | Review Type | Date Requested | Status |
---|---|---|---|
Nick Papior | Needs Fixing | ||
Yann Pouillon | Approve | ||
Review via email: mp+364246@code.launchpad.net |
Commit message
Description of the change
Decouple MPI from fdf.
For an example of usage in Siesta, see the psml-fdf branch.
I did not merge Nick's fix-filedebug branch yet, but most of the changes are no longer
needed as multi-node logging is now done by client programs (we can still advise users
on this).
Note also that I have not yet included the 'broadcast_
in the library (it is in the Siesta branch referenced above). What is the best place to put it?
This merge request is thus work-in-progress, to clean the building system of MPI issues, etc, and
to streamline things. I bumped the version to 0.2.0.
Nick Papior (nickpapior) wrote : | # |
Looks great! I have some minor things, mainly naming convention fdf_* rather than *_fdf_*.
Also some very minor things for consistency, character(len=1) rather than character? (completely optional)
Otherwise approved!
Preview Diff
1 | diff --git a/README.md b/README.md |
2 | index a350e4b..9d49622 100644 |
3 | --- a/README.md |
4 | +++ b/README.md |
5 | @@ -6,14 +6,15 @@ for Electronic Structure Calculations. |
6 | |
7 | LibFDF is the official implementation of the FDF Specifications. |
8 | |
9 | -The LibFDF package, designed by Alberto Garcia and Jose Soler, has been |
10 | -reimplemented by Raul de la Cruz (Barcelona Supercomputer Center) in 2007. |
11 | +The LibFDF package, designed by Alberto Garcia and Jose Soler, was |
12 | +reimplemented by Raul de la Cruz (Barcelona Supercomputer Center) in 2007 |
13 | +to make it memory-side instead of file-oriented. |
14 | + |
15 | +The MPI-enabled library is now being turned into a single-node smaller |
16 | +library with the basic functionality. MPI operation will be delegated |
17 | +to the client code, which will just have to coordinate the broadcast |
18 | +of the internal data structure containing the fdf data. |
19 | |
20 | -The library is memory-side instead of file-oriented, and works better |
21 | -and faster than before. Supports MPI executions (outputs have |
22 | -different names for each node), Cluster mode (at least one node only |
23 | -has access to the input file) and Blocking mode (reading phase is |
24 | -split in several steps for huge execution in shared-filesystems). |
25 | |
26 | We request that the copyright notices in the code be kept if the package |
27 | is eventually used as part of another program. |
28 | diff --git a/configure.ac b/configure.ac |
29 | index 23307c5..a45483c 100644 |
30 | --- a/configure.ac |
31 | +++ b/configure.ac |
32 | @@ -12,7 +12,7 @@ configuring-libfdf.md file.]) |
33 | |
34 | # Init Autoconf |
35 | AC_PREREQ([2.69]) |
36 | -AC_INIT([LibFDF], [0.1.1], |
37 | +AC_INIT([LibFDF], [0.2.1], |
38 | [https://bugs.launchpad.net/libfdf], |
39 | [libfdf], |
40 | [https://launchpad.net/libfdf]) |
41 | diff --git a/doc/fdf.Standard b/doc/fdf.Standard |
42 | index f0cbc72..98703fd 100644 |
43 | --- a/doc/fdf.Standard |
44 | +++ b/doc/fdf.Standard |
45 | @@ -68,12 +68,6 @@ Label < Filename # Comments |
46 | |
47 | By default, all the fdf requests are logged, printing the final value |
48 | extracted (if it is the default value, it is identified as such). |
49 | -In MPI operation, only the master node does the logging, unless the |
50 | -"output-level" is set to a value >= 2. This can be done directly in |
51 | -the fdf file: |
52 | - |
53 | -fdf-output 2 # Turn on logging for all nodes |
54 | -fdf-output 0 # Turn off logging completely |
55 | |
56 | No debugging is done, unless the "debug-level" is set to a value greater |
57 | than zero. The most meaningful way to use this feature is to set the |
58 | @@ -82,14 +76,3 @@ achieved by calling the routine fdf_setdebug with the appropriate level |
59 | and file name *before* the call to fdf_init. This gives full control |
60 | over the behavior. |
61 | |
62 | -Alternatively, if the library is compiled with the pre-processor symbol |
63 | -FDF_DEBUG defined, it will set the debug level to 2 (exhaustive) for all |
64 | -nodes in the system. (It will also set the output level to 2 (all nodes)). |
65 | - |
66 | -If the debugging level (>=2) is only specified in the fdf file itself, |
67 | -the library will provide a print-out of the final data structure generated |
68 | -(in all nodes). |
69 | - |
70 | -fdf-debug 2 |
71 | - |
72 | - |
73 | diff --git a/src/fdf.F90 b/src/fdf.F90 |
74 | index 475ec12..eef0cd7 100644 |
75 | --- a/src/fdf.F90 |
76 | +++ b/src/fdf.F90 |
77 | @@ -26,6 +26,7 @@ |
78 | ! do while(fdf_bline(bfdf, pline)) |
79 | ! (process line 'integers|reals|values|names ...') |
80 | ! enddo |
81 | +! call fdf_bclose(bfdf) |
82 | ! endif |
83 | ! |
84 | ! The subroutine 'fdf_block' returns in 'bfdf' a structure used |
85 | @@ -35,6 +36,9 @@ |
86 | ! line, non-comment line from the block, unless there are no more |
87 | ! lines, in which case it returns .FALSE. and 'pline' is undefined. |
88 | ! |
89 | +! Routine fdf_bclose runs the remaining lines in the block and ensures |
90 | +! the log may be used as input in subsequent entries. |
91 | +! |
92 | ! Routine 'backspace' moves the internal pointer of 'block_fdf' |
93 | ! structure to the previous line returned. |
94 | ! |
95 | @@ -58,37 +62,6 @@ |
96 | ! thread-safe because each thread keeps its relative information |
97 | ! about the search/query that the caller program requests. |
98 | ! |
99 | -! 2) MPI-aware: For MPI executions, FDF library renames output/debugging |
100 | -! or log files to prevent overlaps of these files in a |
101 | -! shared/parallel filesystem. |
102 | -! |
103 | -! This option is enabled with the HAVE_MPI macro, and superseded by |
104 | -! the HAVE_CLUSTER_IO and HAVE_BLOCKING_IO macros. |
105 | -! |
106 | -! 3) Cluster filesystem: It is able to read the input FDF file from a |
107 | -! non-shared filesystem and broadcast the information to the rest |
108 | -! of the nodes in the execution. |
109 | -! |
110 | -! The implementation is as follows: the first node in the MPI rank |
111 | -! that is the owner of the input FDF file, process the file and |
112 | -! send it to the rest of the nodes in the MPI communicator. |
113 | -! |
114 | -! This option is enabled with HAVE_CLUSTER_IO macro. This option cannot be |
115 | -! used if HAVE_BLOCKING_IO macro is enabled. |
116 | -! |
117 | -! 4) Blocking reading: For huge executions (> 1.000 nodes) this option |
118 | -! is useful for shared/parallel filesystems where reading the same |
119 | -! file by several nodes could be a problem (collapsing system). |
120 | -! |
121 | -! The implementation is as follows: reading phase is done in blocking |
122 | -! pattern of size FDF_BLOCKSIZE nodes (configurable through the |
123 | -! FDF_BLOCKSIZE environment variable). |
124 | -! This means that the number of steps needed is STEPS = |
125 | -! #NODES / FDF_BLOCKSIZE. |
126 | -! |
127 | -! This option is enabled with HAVE_BLOCKING_IO macro. This option |
128 | -! cannot be used if HAVE_CLUSTER_IO macro is enabled. |
129 | -! |
130 | ! Alberto Garcia, 1996-2007 |
131 | ! Raul de la Cruz (BSC), September 2007 |
132 | ! |
133 | @@ -122,14 +95,14 @@ MODULE fdf |
134 | |
135 | ! User callable routines in FDF library |
136 | |
137 | -! Start, stop and check parallelism in FDF system |
138 | - public :: fdf_init, fdf_shutdown, fdf_parallel |
139 | +! Start, stop FDF system |
140 | + public :: fdf_init, fdf_shutdown |
141 | |
142 | ! Reading label functions |
143 | public :: fdf_get |
144 | public :: fdf_integer, fdf_single, fdf_double |
145 | public :: fdf_string, fdf_boolean |
146 | - public :: fdf_physical, fdf_isphysical, fdf_convfac |
147 | + public :: fdf_physical, fdf_convfac |
148 | |
149 | ! Lists |
150 | public :: fdf_islist, fdf_islinteger, fdf_islreal |
151 | @@ -139,7 +112,7 @@ MODULE fdf |
152 | public :: fdf_getline |
153 | |
154 | ! Test if label is defined |
155 | - public :: fdf_defined |
156 | + public :: fdf_defined, fdf_isphysical, fdf_isblock |
157 | |
158 | ! Allow to overwrite things in the FDF |
159 | public :: fdf_overwrite, fdf_removelabel, fdf_addline |
160 | @@ -148,13 +121,14 @@ MODULE fdf |
161 | public :: fdf_deprecated, fdf_obsolete |
162 | |
163 | ! %block reading (processing each line) |
164 | - public :: fdf_block, fdf_bline, fdf_bbackspace, fdf_brewind |
165 | + public :: fdf_block, fdf_block_linecount |
166 | + public :: fdf_bline, fdf_bbackspace, fdf_brewind, fdf_bclose |
167 | public :: fdf_bnintegers, fdf_bnreals, fdf_bnvalues, fdf_bnnames, fdf_bntokens |
168 | public :: fdf_bintegers, fdf_breals, fdf_bvalues, fdf_bnames, fdf_btokens |
169 | public :: fdf_bboolean, fdf_bphysical |
170 | public :: fdf_bnlists, fdf_bnilists, fdf_bnrlists, fdf_bnvlists |
171 | public :: fdf_bilists, fdf_brlists, fdf_bvlists |
172 | - |
173 | + |
174 | ! Match, search over blocks, and destroy block structure |
175 | public :: fdf_bmatch, fdf_bsearch, fdf_substring_search |
176 | |
177 | @@ -174,21 +148,8 @@ MODULE fdf |
178 | ! Destroy dynamic list of FDF structure (called in fdf_shutdown) |
179 | private :: fdf_destroy, fdf_destroy_dl |
180 | |
181 | -! MPI init/finalize functions |
182 | -#ifdef HAVE_MPI |
183 | - private :: fdf_mpi_init, fdf_mpi_finalize |
184 | -#endif |
185 | - |
186 | -! Reading functions for HAVE_CLUSTER_IO and HAVE_BLOCKING_IO configuration |
187 | -#ifdef HAVE_CLUSTER_IO |
188 | - private :: setup_fdf_cluster, broadcast_fdf_struct |
189 | -#endif |
190 | -#ifdef HAVE_BLOCKING_IO |
191 | - private :: fdf_readblocking |
192 | -#endif |
193 | - |
194 | ! Debugging functions, level and prints debugging info |
195 | - private :: fdf_printfdf |
196 | + public :: fdf_printfdf |
197 | |
198 | ! Finds a label in the FDF herarchy |
199 | private :: fdf_locate |
200 | @@ -300,37 +261,36 @@ MODULE fdf |
201 | integer(ip), private :: fdf_in(maxdepth) |
202 | integer(ip), private :: fdf_out, fdf_err, fdf_log |
203 | |
204 | -! MPI variables (id, number of MPI tasks) |
205 | -#ifdef HAVE_MPI |
206 | - logical, private :: mpiflag |
207 | - integer(ip), private :: rank, ntasks |
208 | -#endif |
209 | - |
210 | ! Structure for searching inside fdf blocks |
211 | type, public :: block_fdf |
212 | character(len=MAX_LENGTH) :: label |
213 | - type(line_dlist), pointer :: mark |
214 | + type(line_dlist), pointer :: mark => null() |
215 | end type block_fdf |
216 | |
217 | ! Dynamic list for parsed_line structures |
218 | type, public :: line_dlist |
219 | character(len=MAX_LENGTH) :: str |
220 | - type(parsed_line), pointer :: pline |
221 | - type(line_dlist), pointer :: next |
222 | - type(line_dlist), pointer :: prev |
223 | + type(parsed_line), pointer :: pline => null() |
224 | + ! |
225 | + type(line_dlist), pointer :: next => null() |
226 | + type(line_dlist), pointer :: prev => null() |
227 | end type line_dlist |
228 | |
229 | ! FDF data structure (first and last lines) |
230 | type, private :: fdf_file |
231 | integer(ip) :: nlines |
232 | - type(line_dlist), pointer :: first |
233 | - type(line_dlist), pointer :: last |
234 | + type(line_dlist), pointer :: first => null() |
235 | + type(line_dlist), pointer :: last => null() |
236 | end type fdf_file |
237 | |
238 | ! Input FDF file |
239 | - type(fdf_file), pointer, private :: file_in |
240 | + type(fdf_file), private :: file_in |
241 | |
242 | +! Export the following to enable serialization by clients of the library |
243 | |
244 | + public :: fdf_serialize_struct |
245 | + public :: fdf_recreate_struct |
246 | + public :: fdf_set_started |
247 | |
248 | ! Define by default all the others inherit module entities as privated |
249 | ! avoiding redefinitions of entities in several module files with same name |
250 | @@ -339,7 +299,7 @@ MODULE fdf |
251 | private |
252 | |
253 | |
254 | - CONTAINS |
255 | +CONTAINS |
256 | |
257 | ! |
258 | ! Initialization for fdf. |
259 | @@ -365,9 +325,6 @@ MODULE fdf |
260 | THIS_FILE, __LINE__, fdf_err) |
261 | endif |
262 | |
263 | -#ifdef HAVE_MPI |
264 | - call fdf_mpi_init() |
265 | -#endif |
266 | call fdf_initdata() |
267 | |
268 | call io_geterr(fdf_err) |
269 | @@ -377,20 +334,10 @@ MODULE fdf |
270 | fileInput, fileOutput, unitInput ) |
271 | filedebug = trim(fileout) // ".debug" |
272 | |
273 | -#ifdef FDF_DEBUG |
274 | - ! To monitor the parsing and the build-up of the |
275 | - ! fdf data structures in all nodes |
276 | - call fdf_setdebug(2,filedebug) |
277 | - call fdf_setoutput(2,fileout) ! All nodes print output |
278 | -#endif |
279 | - |
280 | -!! call fdf_set_output_file(fileout) |
281 | - |
282 | call fdf_input(filein) |
283 | |
284 | fdf_started = .TRUE. |
285 | |
286 | -#ifndef FDF_DEBUG |
287 | ! Flags within the fdf file itself. |
288 | |
289 | ! At this point only the final fdf data structure will be shown, |
290 | @@ -401,19 +348,14 @@ MODULE fdf |
291 | ! The default is to have output only in the master node |
292 | output_level = fdf_get('fdf-output', 1) |
293 | call fdf_setoutput(output_level,fileout) |
294 | -#endif |
295 | |
296 | - if (fdf_debug2) call fdf_printfdf() |
297 | + if (debug_level >= 2) call fdf_printfdf() |
298 | |
299 | -#ifdef HAVE_MPI |
300 | - call fdf_mpi_finalize() |
301 | -#endif |
302 | !$OMP END SINGLE |
303 | RETURN |
304 | !------------------------------------------------------------------------- END |
305 | END SUBROUTINE fdf_init |
306 | |
307 | - |
308 | SUBROUTINE set_file_names( fileIn, fileOut, & |
309 | optFileIn, optFileOut, unitIn ) |
310 | ! If present, copies input arguments optFileIn/Out to fileIn/Out. |
311 | @@ -433,14 +375,10 @@ MODULE fdf |
312 | |
313 | integer:: count, ierr, iostat, iu, iuIn |
314 | logical:: named, opened |
315 | - character(len=300) line |
316 | + character(len=MAX_LENGTH*2) line |
317 | character(len=maxFileNameLength) fileName |
318 | |
319 | !------------------------------------------------------------------------- BEGIN |
320 | -#ifdef HAVE_MPI |
321 | - if (rank==0) then |
322 | -#endif |
323 | - |
324 | ! Find a job-specific number |
325 | call system_clock( count ) |
326 | count = mod(count,100000) |
327 | @@ -507,83 +445,11 @@ MODULE fdf |
328 | endif ! (fileName=='stdin') |
329 | |
330 | endif ! (present(optFileIn)) |
331 | - |
332 | -#ifdef HAVE_MPI |
333 | - endif ! (rank==0) |
334 | -#endif |
335 | !--------------------------------------------------------------------------- END |
336 | END SUBROUTINE set_file_names |
337 | |
338 | ! |
339 | -! Initialize MPI subsystem if the application calling/using FDF |
340 | -! library is not running with MPI enabled. |
341 | -! |
342 | -#ifdef HAVE_MPI |
343 | - SUBROUTINE fdf_mpi_init() |
344 | - use mpi_siesta |
345 | - implicit none |
346 | -!--------------------------------------------------------------- Local Variables |
347 | - integer(ip) :: ierr |
348 | - |
349 | -!------------------------------------------------------------------------- BEGIN |
350 | - call MPI_Initialized(mpiflag, ierr) |
351 | - if (ierr .ne. MPI_SUCCESS) then |
352 | - call die('FDF module: fdf_mpi_init', 'Error initializing MPI system.' // & |
353 | - 'Terminating.', THIS_FILE, __LINE__, fdf_err, rc=ierr) |
354 | - endif |
355 | - |
356 | - if (.not. mpiflag) then |
357 | - call MPI_Init(ierr) |
358 | - if (ierr .ne. MPI_SUCCESS) then |
359 | - call die('FDF module: fdf_mpi_init', 'Error initializing MPI system.' // & |
360 | - 'Terminating.', THIS_FILE, __LINE__, fdf_err, rc=ierr) |
361 | - endif |
362 | - endif |
363 | - |
364 | - call MPI_Comm_Rank(MPI_COMM_WORLD, rank, ierr) |
365 | - if (ierr .ne. MPI_SUCCESS) then |
366 | - call die('FDF module: fdf_mpi_init', 'Error getting MPI comm rank.' // & |
367 | - 'Terminating.', THIS_FILE, __LINE__, fdf_err, rc= ierr) |
368 | - endif |
369 | - |
370 | - call MPI_Comm_Size(MPI_COMM_WORLD, ntasks, ierr) |
371 | - if (ierr .ne. MPI_SUCCESS) then |
372 | - call die('FDF module: fdf_mpi_init', 'Error getting MPI comm size.' // & |
373 | - 'Terminating.', THIS_FILE, __LINE__, fdf_err, rc=ierr) |
374 | - endif |
375 | - RETURN |
376 | -!--------------------------------------------------------------------------- END |
377 | - END SUBROUTINE fdf_mpi_init |
378 | -#endif |
379 | - |
380 | -! |
381 | -! Finalize MPI subsystem if the application calling/using FDF |
382 | -! library is not running with MPI enabled. |
383 | -! |
384 | -#ifdef HAVE_MPI |
385 | - SUBROUTINE fdf_mpi_finalize() |
386 | - use mpi_siesta |
387 | - implicit none |
388 | -!--------------------------------------------------------------- Local Variables |
389 | - integer(ip) :: ierr |
390 | - |
391 | -!------------------------------------------------------------------------- BEGIN |
392 | - if (.not. mpiflag) then |
393 | - call MPI_Finalize(ierr) |
394 | - if (ierr .ne. MPI_SUCCESS) then |
395 | - call die('FDF module: fdf_mpi_finalize', 'Error finalizing MPI system.' // & |
396 | - 'Terminating.', THIS_FILE, __LINE__, fdf_err, rc=ierr) |
397 | - endif |
398 | - endif |
399 | - |
400 | - RETURN |
401 | -!--------------------------------------------------------------------------- END |
402 | - END SUBROUTINE fdf_mpi_finalize |
403 | -#endif |
404 | - |
405 | -! |
406 | -! Reads the input file depending on the configuration of the system: |
407 | -! Shared Filesystem, Cluster Filesystem or Blocking input access |
408 | +! Reads the input file |
409 | ! |
410 | SUBROUTINE fdf_input(filein) |
411 | implicit none |
412 | @@ -591,14 +457,8 @@ MODULE fdf |
413 | character(*) :: filein |
414 | |
415 | !------------------------------------------------------------------------- BEGIN |
416 | -#ifdef HAVE_CLUSTER_IO |
417 | -!! call fdf_readcluster(filein) |
418 | - call setup_fdf_cluster(filein) |
419 | -#elif defined(HAVE_BLOCKING_IO) |
420 | - call fdf_readblocking(filein) |
421 | -#else |
422 | + |
423 | call fdf_read(filein) |
424 | -#endif |
425 | |
426 | if (fdf_output) write(fdf_out,'(a,a,a,i3)') '#FDF module: Opened ', filein, & |
427 | ' for input. Unit:', fdf_in(1) |
428 | @@ -607,235 +467,6 @@ MODULE fdf |
429 | !--------------------------------------------------------------------------- END |
430 | END SUBROUTINE fdf_input |
431 | |
432 | -! |
433 | -! Reading code for Cluster Architecture, where the filesystem |
434 | -! is not shared along all the nodes in the system. |
435 | -! |
436 | -#ifdef HAVE_CLUSTER_IO |
437 | - SUBROUTINE fdf_readcluster(filein) |
438 | - use mpi_siesta |
439 | - implicit none |
440 | -!--------------------------------------------------------------- Input Variables |
441 | - character(*) :: filein |
442 | - |
443 | -!--------------------------------------------------------------- Local Variables |
444 | - character(80) :: msg |
445 | - character(256) :: fileinTmp |
446 | - integer(ip) :: ierr, texist_send, texist_recv |
447 | -#ifdef SOPHISTICATED_SEARCH |
448 | - logical :: file_exist |
449 | -#endif |
450 | -!------------------------------------------------------------------------- BEGIN |
451 | -! Tests if the running node has the input file: |
452 | -! If found: texist_send = rank |
453 | -! Else : texist_send = error_code (ntasks + 1) |
454 | - |
455 | -#ifdef SOPHISTICATED_SEARCH |
456 | - INQUIRE(file=filein, exist=file_exist) |
457 | - if (file_exist) then |
458 | - texist_send = rank |
459 | - else |
460 | - texist_send = ntasks + 1 |
461 | - endif |
462 | - |
463 | - call MPI_AllReduce(texist_send, texist_recv, 1, MPI_INTEGER, & |
464 | - MPI_MIN, MPI_COMM_WORLD, ierr) |
465 | - if (ierr .ne. MPI_SUCCESS) then |
466 | - call die('FDF module: fdf_readcluster', 'Error in MPI_AllReduce (task_exist).' // & |
467 | - 'Terminating.', THIS_FILE, __LINE__, fdf_err, rc=ierr) |
468 | - endif |
469 | -#else |
470 | -! |
471 | -! Simplify: Assume node 0 has the file |
472 | -! |
473 | - texist_recv = 0 |
474 | -#endif |
475 | - |
476 | -! The node owner of the input file send the data to the other ones |
477 | -! if none node has an input file with such name abort the application |
478 | - if (texist_recv .eq. ntasks + 1) then |
479 | - write(msg,*) 'No node found in the cluster with ', & |
480 | - ' input file ', filein,'. Terminating.' |
481 | - call die('FDF module: fdf_readcluster', msg, & |
482 | - THIS_FILE, __LINE__, fdf_err, rc=ierr) |
483 | - else |
484 | - if (texist_recv .eq. rank) then |
485 | - call fdf_read(filein) |
486 | - call fdf_sendInput() |
487 | -!!Debug call MPI_Barrier( MPI_COMM_WORLD, ierr ) |
488 | - if (fdf_output) write(fdf_out,*) '#FDF module: Node', rank, 'reading/sending', & |
489 | - ' input file ', filein |
490 | - else |
491 | - call fdf_recvInput(texist_recv, filein, fileinTmp) |
492 | -!!Debug call MPI_Barrier( MPI_COMM_WORLD, ierr ) |
493 | - call fdf_read(fileinTmp) |
494 | - if (fdf_output) write(fdf_out,*) '#FDF module: Node', rank, 'receiving input', & |
495 | - ' file from', texist_recv, 'to ', TRIM(fileinTmp) |
496 | - endif |
497 | - endif |
498 | - RETURN |
499 | -!--------------------------------------------------------------------------- END |
500 | - END SUBROUTINE fdf_readcluster |
501 | -! |
502 | -! |
503 | - SUBROUTINE setup_fdf_cluster(filein) |
504 | -! |
505 | -! A more efficient alternative to fdf_sendInput/fdf_recvInput |
506 | -! that avoids the creation of scratch files by non-root nodes. |
507 | -! |
508 | -! Alberto Garcia, April 2011 |
509 | - |
510 | - use mpi_siesta |
511 | - implicit none |
512 | -!--------------------------------------------------------------- Input Variables |
513 | - character(*), intent(in) :: filein ! File name |
514 | - |
515 | -!--------------------------------------------------------------- Local Variables |
516 | - character(80) :: msg |
517 | - character(256) :: fileinTmp |
518 | - integer(ip) :: ierr, texist_send, reading_node |
519 | -#ifdef SOPHISTICATED_SEARCH |
520 | - logical :: file_exist |
521 | -#endif |
522 | -!------------------------------------------------------------------------- BEGIN |
523 | -! Tests if the running node has the input file: |
524 | -! If found: texist_send = rank |
525 | -! Else : texist_send = error_code (ntasks + 1) |
526 | - |
527 | -#ifdef SOPHISTICATED_SEARCH |
528 | - INQUIRE(file=filein, exist=file_exist) |
529 | - if (file_exist) then |
530 | - texist_send = rank |
531 | - else |
532 | - texist_send = ntasks + 1 |
533 | - endif |
534 | - |
535 | - call MPI_AllReduce(texist_send, reading_node, 1, MPI_INTEGER, & |
536 | - MPI_MIN, MPI_COMM_WORLD, ierr) |
537 | - if (ierr .ne. MPI_SUCCESS) then |
538 | - call die('FDF module: fdf_readcluster', 'Error in MPI_AllReduce (task_exist).' // & |
539 | - 'Terminating.', THIS_FILE, __LINE__, fdf_err, rc=ierr) |
540 | - endif |
541 | -#else |
542 | -! |
543 | -! Simplify: Assume node 0 has the file |
544 | -! |
545 | - reading_node = 0 |
546 | -#endif |
547 | - |
548 | -! If no node has an input file with such name abort the application |
549 | - if (reading_node .eq. ntasks + 1) then |
550 | - write(msg,*) 'No node found in the cluster with ', & |
551 | - ' input file ', filein,'. Terminating.' |
552 | - call die('FDF module: fdf_readcluster', msg, & |
553 | - THIS_FILE, __LINE__, fdf_err, rc=ierr) |
554 | - endif |
555 | - |
556 | -! The root node reads and digests the input file |
557 | -! and sends the serialized fdf_file data structure to the other nodes |
558 | - |
559 | - if (rank == reading_node) then |
560 | - call fdf_read(filein) |
561 | - endif |
562 | - call broadcast_fdf_struct(reading_node) |
563 | - |
564 | - RETURN |
565 | -!--------------------------------------------------------------------------- END |
566 | - END SUBROUTINE setup_fdf_cluster |
567 | - |
568 | -! |
569 | -! Broadcast complete fdf structure |
570 | -! |
571 | - SUBROUTINE broadcast_fdf_struct(reading_node) |
572 | - use mpi_siesta |
573 | - implicit none |
574 | - |
575 | - integer, intent(in) :: reading_node ! Node which contains the struct |
576 | - |
577 | -!--------------------------------------------------------------- Local Variables |
578 | - character, pointer :: bufferFDF(:) => null() |
579 | - integer(ip) :: i, j, k, ierr, nlines |
580 | - |
581 | -!------------------------------------------------------------------------- BEGIN |
582 | - |
583 | - if (rank == reading_node) then |
584 | - nlines = file_in%nlines |
585 | - endif |
586 | - |
587 | - call MPI_Bcast(nlines, 1, & |
588 | - MPI_INTEGER, reading_node, MPI_COMM_WORLD, ierr) |
589 | - if (ierr .ne. MPI_SUCCESS) then |
590 | - call die('FDF module: broadcast_fdf', 'Error Broadcasting nlines.' // & |
591 | - 'Terminating.', THIS_FILE, __LINE__, fdf_err, rc=ierr) |
592 | - endif |
593 | - |
594 | - ALLOCATE(bufferFDF(nlines*SERIALIZED_LENGTH), stat=ierr) |
595 | - if (ierr .ne. 0) then |
596 | - call die('FDF module: broadcast_fdf', 'Error allocating bufferFDF', & |
597 | - THIS_FILE, __LINE__, fdf_err, rc=ierr) |
598 | - endif |
599 | - |
600 | - if (rank == reading_node) then |
601 | - call serialize_fdf_struct(bufferFDF) |
602 | - endif |
603 | - |
604 | - call MPI_Bcast(bufferFDF, size(bufferFDF), & |
605 | - MPI_CHARACTER, reading_node, MPI_COMM_WORLD, ierr) |
606 | - if (ierr .ne. MPI_SUCCESS) then |
607 | - call die('FDF module: broadcast_fdf', 'Error Broadcasting bufferFDF.' // & |
608 | - 'Terminating.', THIS_FILE, __LINE__, fdf_err, rc=ierr) |
609 | - endif |
610 | - if (rank /= reading_node) then |
611 | - call recreate_fdf_struct(nlines,bufferFDF) |
612 | - endif |
613 | - |
614 | - DEALLOCATE(bufferFDF) |
615 | - |
616 | - RETURN |
617 | -!--------------------------------------------------------------------------- END |
618 | - END SUBROUTINE broadcast_fdf_struct |
619 | -#endif |
620 | -! |
621 | -! Reading code for Blocking read. The reading of the input file is |
622 | -! splitted in several steps of size BLOCKSIZE. The number of steps |
623 | -! needed to read the file in all the nodes depends on the total nodes, |
624 | -! being ntasks/HAVE_BLOCKING_IOSIZE. |
625 | -! |
626 | -! With this method we can avoid a colapse of a global parallel filesystem |
627 | -! if the execution of the application runs over a huge amount of nodes. |
628 | -! |
629 | -#ifdef HAVE_BLOCKING_IO |
630 | - SUBROUTINE fdf_readblocking(filein) |
631 | - use mpi_siesta |
632 | - implicit none |
633 | -!--------------------------------------------------------------- Input Variables |
634 | - character(*) :: filein |
635 | - |
636 | -!--------------------------------------------------------------- Local Variables |
637 | - integer(ip) :: i, ierr |
638 | - |
639 | -!------------------------------------------------------------------------- BEGIN |
640 | - |
641 | - do i= 0, ntasks-1, FDF_BLOCKSIZE |
642 | - if ((rank .ge. i) .and. (rank .le. i+FDF_BLOCKSIZE-1)) then |
643 | - call fdf_read(filein) |
644 | - if (fdf_output) write(fdf_out,*) '#FDF module: Task', rank, 'reading input', & |
645 | - ' file ', filein, ' in step', (i/FDF_BLOCKSIZE)+1 |
646 | - endif |
647 | - |
648 | - call MPI_Barrier(MPI_COMM_WORLD, ierr) |
649 | - if (ierr .ne. MPI_SUCCESS) then |
650 | - call die('FDF module: fdf_readblocking', 'Error in MPI_Barrier (fdf_read).' // & |
651 | - 'Terminating.', THIS_FILE, __LINE__, fdf_err, rc=ierr) |
652 | - endif |
653 | - enddo |
654 | - |
655 | - RETURN |
656 | -!--------------------------------------------------------------------------- END |
657 | - END SUBROUTINE fdf_readblocking |
658 | -#endif |
659 | - |
660 | -! |
661 | ! Read an input file (and include files) and builds memory |
662 | ! structure that will contain the data and will help in searching |
663 | ! |
664 | @@ -1431,12 +1062,6 @@ MODULE fdf |
665 | !------------------------------------------------------------------------- BEGIN |
666 | ndepth = 0 |
667 | |
668 | - ALLOCATE(file_in, stat=ierr) |
669 | - if (ierr .ne. 0) then |
670 | - call die('FDF module: fdf_initdata', 'Error allocating file_in', & |
671 | - THIS_FILE, __LINE__, fdf_err, rc=ierr) |
672 | - endif |
673 | - |
674 | file_in%nlines = 0 |
675 | NULLIFY(file_in%first) |
676 | NULLIFY(file_in%last) |
677 | @@ -1498,7 +1123,7 @@ MODULE fdf |
678 | endif |
679 | |
680 | ! To circumvent the first/last line in the fdf-file |
681 | - ! we have to check for the existance of the |
682 | + ! we have to check for the existence of the |
683 | ! first/last mark being the one removed. |
684 | ! That special case *must* correct the first/last |
685 | ! tokens. |
686 | @@ -1554,7 +1179,7 @@ MODULE fdf |
687 | type(parsed_line), pointer :: pline |
688 | |
689 | !--------------------------------------------------------------- Local Variables |
690 | - integer :: i, ierr |
691 | + integer(ip) :: i, ierr |
692 | type(line_dlist), pointer :: mark |
693 | |
694 | !------------------------------------------------------------------------- BEGIN |
695 | @@ -1670,45 +1295,8 @@ MODULE fdf |
696 | !----------------------------------------------------- BEGIN |
697 | call io_assign(fdf_out) |
698 | |
699 | -#ifdef HAVE_MPI |
700 | - if (rank /= 0) then |
701 | - fileouttmp = "/tmp/" // trim(fileout) // "." // i2s(rank) |
702 | - else |
703 | - fileouttmp = fileout |
704 | - endif |
705 | -#else |
706 | - fileouttmp = fileout |
707 | -#endif |
708 | - |
709 | -#ifdef FDF_DEBUG |
710 | - ! |
711 | - ! If debugging, all the nodes use named log files |
712 | - ! |
713 | - open( unit=fdf_out, file=TRIM(fileouttmp), form='formatted', & |
714 | + open( unit=fdf_out, file=TRIM(fileout), form='formatted', & |
715 | access='sequential', status='replace' ) |
716 | -#else |
717 | - ! |
718 | - ! Only the master node opens a named log file in the current dir. |
719 | - ! Non-master nodes use a scratch file. |
720 | - ! These log files tend to be quite small, so there |
721 | - ! should not be problems such as filling up filesystems. |
722 | - ! ... your mileage might vary. This is a grey area. |
723 | - ! Some compilers allow the user to specify where scratch |
724 | - ! files go. |
725 | - ! |
726 | -#ifdef HAVE_MPI |
727 | - if (rank /= 0) then |
728 | - open( unit=fdf_out, form='formatted', & |
729 | - access='sequential', status='scratch' ) |
730 | - else |
731 | - open( unit=fdf_out, file=TRIM(fileouttmp), form='formatted', & |
732 | - access='sequential', status='replace' ) |
733 | - end if |
734 | -#else |
735 | - open( unit=fdf_out, file=TRIM(fileouttmp), form='formatted', & |
736 | - access='sequential', status='replace' ) |
737 | -#endif |
738 | -#endif |
739 | |
740 | RETURN |
741 | !--------------------------------------------------------------------------- END |
742 | @@ -1739,11 +1327,10 @@ MODULE fdf |
743 | SUBROUTINE fdf_destroy(fdfp) |
744 | implicit none |
745 | !-------------------------------------------------------------- Output Variables |
746 | - type(fdf_file), pointer :: fdfp |
747 | + type(fdf_file) :: fdfp |
748 | |
749 | !------------------------------------------------------------------------- BEGIN |
750 | if (ASSOCIATED(fdfp%first)) call fdf_destroy_dl(fdfp%first) |
751 | - DEALLOCATE(fdfp) |
752 | |
753 | RETURN |
754 | !--------------------------------------------------------------------------- END |
755 | @@ -1757,7 +1344,10 @@ MODULE fdf |
756 | !-------------------------------------------------------------- Output Variables |
757 | type(line_dlist), pointer :: dlp |
758 | |
759 | -!------------------------------------------------------------------------- BEGIN |
760 | + !! Use for tail recursion later: type(line_dlist), pointer :: pnext |
761 | + |
762 | + !------------------------------------------------------------------------- BEGIN |
763 | + ! This is NOT tail-recursive!! |
764 | if (ASSOCIATED(dlp%next)) call fdf_destroy_dl(dlp%next) |
765 | call destroy(dlp%pline) |
766 | DEALLOCATE(dlp) |
767 | @@ -1765,27 +1355,6 @@ MODULE fdf |
768 | RETURN |
769 | !--------------------------------------------------------------------------- END |
770 | END SUBROUTINE fdf_destroy_dl |
771 | - |
772 | -! |
773 | -! Tells when FDF library has been compiled with parallel |
774 | -! execution support. Depends on HAVE_MPI macro. |
775 | -! |
776 | - FUNCTION fdf_parallel() |
777 | - implicit none |
778 | -!-------------------------------------------------------------- Output Variables |
779 | - logical :: fdf_parallel |
780 | - |
781 | -!------------------------------------------------------------------------- BEGIN |
782 | -#ifdef HAVE_MPI |
783 | - fdf_parallel = .TRUE. |
784 | -#else |
785 | - fdf_parallel = .FALSE. |
786 | -#endif |
787 | - |
788 | - RETURN |
789 | -!--------------------------------------------------------------------------- END |
790 | - END FUNCTION fdf_parallel |
791 | - |
792 | ! |
793 | ! Read a line of the 'ndepth' input file, returning .TRUE. if |
794 | ! there are more lines to read from input file, .FALSE. otherwise. |
795 | @@ -2093,6 +1662,7 @@ MODULE fdf |
796 | call integerlists(mark%pline,1,ni,list) |
797 | end if |
798 | |
799 | + ! find a way to write out the list anyway |
800 | if (fdf_output) write(fdf_out,'(a,5x,i10)') label, lni |
801 | else |
802 | write(msg,*) 'no list value for ', label |
803 | @@ -2813,6 +2383,68 @@ MODULE fdf |
804 | END FUNCTION fdf_locate |
805 | |
806 | ! |
807 | +! Returns true or false whether or not the label 'label' is |
808 | +! a block. |
809 | +! I.e. it returns true if the line has the form bl, if not found, or not bl |
810 | +! it returns false. |
811 | +! |
812 | + FUNCTION fdf_isblock(label) |
813 | + implicit none |
814 | +!--------------------------------------------------------------- Input Variables |
815 | + character(*) :: label |
816 | + |
817 | +!-------------------------------------------------------------- Output Variables |
818 | + logical :: fdf_isblock |
819 | + |
820 | +!--------------------------------------------------------------- Local Variables |
821 | + type(line_dlist), pointer :: mark |
822 | + character(80) :: strlabel |
823 | + |
824 | +!------------------------------------------------------------------------- BEGIN |
825 | +! Prevents using FDF routines without initialize |
826 | + if (.not. fdf_started) then |
827 | + call die('FDF module: fdf_isblock', 'FDF subsystem not initialized', & |
828 | + THIS_FILE, __LINE__, fdf_err) |
829 | + endif |
830 | + |
831 | + fdf_isblock = .false. |
832 | + |
833 | + mark => file_in%first |
834 | + do while ( associated(mark) ) |
835 | + |
836 | +!!$ if ( match(mark%pline, 'l') ) then |
837 | +!!$ strlabel = labels(mark%pline) |
838 | +!!$ |
839 | +!!$ if ( labeleq(strlabel, label, fdf_log) ) then |
840 | +!!$ ! fdf has first-encounter acceptance. |
841 | +!!$ ! I.e. for an input |
842 | +!!$ ! Label_Name 1 |
843 | +!!$ ! %block Label_Name |
844 | +!!$ ! 1 |
845 | +!!$ ! %endblock Label_Name |
846 | +!!$ ! the former will be accepted first. |
847 | +!!$ exit |
848 | +!!$ end if |
849 | + |
850 | + if ( match(mark%pline, 'bl') ) then |
851 | + strlabel = blocks(mark%pline) |
852 | + |
853 | + if ( labeleq(strlabel, label, fdf_log) ) then |
854 | + fdf_isblock = .true. |
855 | + exit |
856 | + end if |
857 | + end if |
858 | + |
859 | + mark => mark%next |
860 | + end do |
861 | + |
862 | + if (fdf_output) write(fdf_out,'(a,5x,l10)') "#:block? " // label, fdf_isblock |
863 | + |
864 | + RETURN |
865 | +!--------------------------------------------------------------------------- END |
866 | + END FUNCTION fdf_isblock |
867 | + |
868 | +! |
869 | ! Searches for block label in the fdf hierarchy. If it appears returns |
870 | ! .TRUE. and leaves block mark pointer positioned at the first line. |
871 | ! Otherwise, it returns .FALSE. and block mark points to NULL. |
872 | @@ -2965,7 +2597,7 @@ MODULE fdf |
873 | if (labeleq(strlabel, bfdf%label, fdf_log)) then |
874 | bfdf%mark => bfdf%mark%prev |
875 | |
876 | - if (fdf_output) write(fdf_out,'(1x,a)') "(Backspace to) " // "|" // & |
877 | + if (fdf_output) write(fdf_out,'(1x,a)') "#:(Backspace to) " // "|" // & |
878 | TRIM(bfdf%mark%str) // "|" |
879 | endif |
880 | |
881 | @@ -2976,14 +2608,14 @@ MODULE fdf |
882 | |
883 | if (labeleq(strlabel, bfdf%label, fdf_log)) then |
884 | fdf_bbackspace = .FALSE. |
885 | - if (fdf_output) write(fdf_out,'(1x,a)') "(Cannot backspace) " // "|" // & |
886 | + if (fdf_output) write(fdf_out,'(1x,a)') "#:(Cannot backspace) " // "|" // & |
887 | TRIM(bfdf%mark%str) // "|" |
888 | endif |
889 | |
890 | else |
891 | |
892 | bfdf%mark => bfdf%mark%prev |
893 | - if (fdf_output) write(fdf_out,'(1x,a)') "(Backspace to) " // "|" // & |
894 | + if (fdf_output) write(fdf_out,'(1x,a)') "#:(Backspace to) " // "|" // & |
895 | TRIM(bfdf%mark%str) // "|" |
896 | endif |
897 | |
898 | @@ -3028,6 +2660,108 @@ MODULE fdf |
899 | END SUBROUTINE fdf_brewind |
900 | |
901 | ! |
902 | +! Closes the opened block by looping the remaining lines of the working line. |
903 | +! This is only needed to complete the fdf-*.log files output for direct |
904 | +! usage later. It does nothing internally. |
905 | +! |
906 | + SUBROUTINE fdf_bclose(bfdf) |
907 | + implicit none |
908 | +!-------------------------------------------------------------- Output Variables |
909 | + type(block_fdf) :: bfdf |
910 | + |
911 | +!--------------------------------------------------------------- Local Variables |
912 | + type(parsed_line), pointer :: pline |
913 | + integer(ip) :: i |
914 | + character(80) :: msg |
915 | + |
916 | +!------------------------------------------------------------------------- BEGIN |
917 | +! Prevents using FDF routines without initialize |
918 | + if (.not. fdf_started) then |
919 | + call die('FDF module: fdf_bclose', 'FDF subsystem not initialized', & |
920 | + THIS_FILE, __LINE__, fdf_err) |
921 | + endif |
922 | + |
923 | + ! Quick return (no need for errors) |
924 | + if ( .not. associated(bfdf%mark) ) return |
925 | + |
926 | + ! This should hopefully discourage compilers to optimize the loop away... |
927 | + i = 0 |
928 | + do while ( fdf_bline(bfdf, pline) ) |
929 | + i = i + fdf_bnvalues(pline) |
930 | + end do |
931 | + write(msg,'(a,i10)') 'Block ', i |
932 | + |
933 | + RETURN |
934 | +!--------------------------------------------------------------------------- END |
935 | + END SUBROUTINE fdf_bclose |
936 | + |
937 | + |
938 | +! |
939 | +! Count number of lines with an optional specification. |
940 | +! I.e. this will read through the block and return the number of lines in the |
941 | +! block which matches the morphology (morph) |
942 | +! This may be used to easily digest number of non-empty lines in the block. |
943 | +! Note that a match on the morphology only compares the number of ID's in |
944 | +! morph. I.e. a line with 'vvvil' will match 'vvvi'. |
945 | +! |
946 | + FUNCTION fdf_block_linecount(label, morph) |
947 | + implicit none |
948 | +!--------------------------------------------------------------- Input Variables |
949 | + character(len=*) :: label |
950 | + character(len=*), optional :: morph |
951 | +!-------------------------------------------------------------- Output Variables |
952 | + integer(ip) :: fdf_block_linecount |
953 | + |
954 | +!--------------------------------------------------------------- Local Variables |
955 | + type(block_fdf) :: bfdf |
956 | + type(parsed_line), pointer :: pline |
957 | + logical :: orig_fdf_output |
958 | + |
959 | +!------------------------------------------------------------------------- BEGIN |
960 | +! Prevents using FDF routines without initialize |
961 | + if (.not. fdf_started) then |
962 | + call die('FDF module: fdf_block_linecount', 'FDF subsystem not initialized', & |
963 | + THIS_FILE, __LINE__, fdf_err) |
964 | + endif |
965 | + |
966 | + ! Store the fdf_output variable (suppress writing to log) |
967 | + orig_fdf_output = fdf_output |
968 | + fdf_output = .false. |
969 | + |
970 | + ! Find the block and search for morhp |
971 | + fdf_block_linecount = 0 |
972 | + if ( fdf_block(label, bfdf) ) then |
973 | + |
974 | + do while ( fdf_bline(bfdf, pline) ) |
975 | + if ( present(morph) ) then |
976 | + if ( fdf_bmatch(pline, morph) ) then |
977 | + fdf_block_linecount = fdf_block_linecount + 1 |
978 | + end if |
979 | + else |
980 | + fdf_block_linecount = fdf_block_linecount + 1 |
981 | + end if |
982 | + end do |
983 | + |
984 | + end if |
985 | + |
986 | + ! Restore output |
987 | + fdf_output = orig_fdf_output |
988 | + |
989 | + if ( fdf_output ) then |
990 | + if ( present(morph) ) then |
991 | + write(fdf_out,'(3a,3x,i0)') "#:block-line-count? ", & |
992 | + trim(label), ' ('//trim(morph)//')', fdf_block_linecount |
993 | + else |
994 | + write(fdf_out,'(2a,3x,i0)') "#:block-line-count? ", & |
995 | + trim(label), fdf_block_linecount |
996 | + end if |
997 | + end if |
998 | + |
999 | + RETURN |
1000 | +!--------------------------------------------------------------------------- END |
1001 | + END FUNCTION fdf_block_linecount |
1002 | + |
1003 | +! |
1004 | ! Check if label is defined |
1005 | ! |
1006 | logical FUNCTION fdf_defined(label) |
1007 | @@ -3037,18 +2771,16 @@ MODULE fdf |
1008 | |
1009 | !--------------------------------------------------------------- Local Variables |
1010 | type(line_dlist), pointer :: mark |
1011 | - type(block_fdf) :: bfdf |
1012 | - |
1013 | |
1014 | !--------------------------------------------------------------------- BEGIN |
1015 | ! First, check whether a single label exists: |
1016 | fdf_defined = fdf_locate(label, mark) |
1017 | if (.not. fdf_defined) then |
1018 | ! Check whether there is a block with that label |
1019 | - fdf_defined = fdf_block(label,bfdf) |
1020 | + fdf_defined = fdf_isblock(label) |
1021 | endif |
1022 | - if (fdf_defined) then |
1023 | - if (fdf_output) write(fdf_out,'(a)') label |
1024 | + if ( fdf_output ) then |
1025 | + write(fdf_out,'(a,5x,l10)') '#:defined? ' // label, fdf_defined |
1026 | endif |
1027 | |
1028 | RETURN |
1029 | @@ -3077,30 +2809,11 @@ MODULE fdf |
1030 | endif |
1031 | else |
1032 | if (.not. fdf_output) then |
1033 | -#ifdef HAVE_MPI |
1034 | - if (rank /= 0) then |
1035 | - if (level .ge. 2) then |
1036 | - fileout = "/tmp/" // trim(fileout) // "." // i2s(rank) |
1037 | - call io_assign(fdf_out) |
1038 | - open(fdf_out, file=fileout, form='formatted', & |
1039 | - status='unknown') |
1040 | - REWIND(fdf_out) |
1041 | - fdf_output = .TRUE. |
1042 | - endif |
1043 | - else |
1044 | - call io_assign(fdf_out) |
1045 | - open(fdf_out, file=fileout, form='formatted', & |
1046 | - status='unknown') |
1047 | - REWIND(fdf_out) |
1048 | - fdf_output = .TRUE. |
1049 | - endif |
1050 | -#else |
1051 | call io_assign(fdf_out) |
1052 | open(fdf_out, file=fileout, form='formatted', & |
1053 | status='unknown') |
1054 | REWIND(fdf_out) |
1055 | fdf_output = .TRUE. |
1056 | -#endif |
1057 | endif |
1058 | endif |
1059 | !----------------------------------------------------------------------- END |
1060 | @@ -3127,11 +2840,6 @@ MODULE fdf |
1061 | else |
1062 | if (.not. fdf_debug) then |
1063 | call io_assign(fdf_log) |
1064 | -#ifdef HAVE_MPI |
1065 | - if (rank /= 0) then |
1066 | - filedebug = "/tmp/" // trim(filedebug) // "." // i2s(rank) |
1067 | - endif |
1068 | -#endif |
1069 | open(fdf_log, file=filedebug, form='formatted', & |
1070 | status='unknown') |
1071 | REWIND(fdf_log) |
1072 | @@ -3162,13 +2870,13 @@ MODULE fdf |
1073 | |
1074 | !------------------------------------------------------------------------- BEGIN |
1075 | if ( fdf_defined(label) ) then |
1076 | - if (fdf_output) write(fdf_out,'(a)') "**Warning: FDF symbol '"//trim(label)// & |
1077 | + if (fdf_output) write(fdf_out,'(a)') "#**Warning: FDF symbol '"//trim(label)// & |
1078 | "' is deprecated." |
1079 | if ( fdf_defined(newlabel) ) then |
1080 | - if (fdf_output) write(fdf_out,'(a)') " FDF symbol '"//trim(newlabel)// & |
1081 | + if (fdf_output) write(fdf_out,'(a)') "# FDF symbol '"//trim(newlabel)// & |
1082 | "' will be used instead." |
1083 | else |
1084 | - if (fdf_output) write(fdf_out,'(a)') " FDF symbol '"//trim(newlabel)// & |
1085 | + if (fdf_output) write(fdf_out,'(a)') "# FDF symbol '"//trim(newlabel)// & |
1086 | "' replaces '"//trim(label)//"'." |
1087 | end if |
1088 | end if |
1089 | @@ -3186,21 +2894,28 @@ MODULE fdf |
1090 | |
1091 | !------------------------------------------------------------------------- BEGIN |
1092 | if ( fdf_defined(label) ) then |
1093 | - if (fdf_output) write(fdf_out,'(a)') "**Warning: FDF symbol '"//trim(label)// & |
1094 | + if (fdf_output) write(fdf_out,'(a)') "#**Warning: FDF symbol '"//trim(label)// & |
1095 | "' is obsolete." |
1096 | end if |
1097 | |
1098 | !--------------------------------------------------------------------------- END |
1099 | end subroutine fdf_obsolete |
1100 | |
1101 | - |
1102 | - subroutine serialize_fdf_struct(buffer) |
1103 | - character, pointer :: buffer(:) |
1104 | +!===================== Serialization utilities for clients |
1105 | + |
1106 | + subroutine fdf_serialize_struct(buffer) |
1107 | + character(len=1), intent(inout), allocatable :: buffer(:) |
1108 | |
1109 | character(len=SERIALIZED_LENGTH) bufline |
1110 | type(line_dlist), pointer :: mark |
1111 | - integer :: i, length, init, final |
1112 | + integer(ip) :: i, length, init, final |
1113 | + |
1114 | + integer :: nchars ! total size of serialized content |
1115 | |
1116 | + if (allocated(buffer)) deallocate(buffer) |
1117 | + nchars = file_in%nlines * SERIALIZED_LENGTH |
1118 | + allocate(buffer(nchars)) |
1119 | + |
1120 | mark => file_in%first |
1121 | do i= 1, file_in%nlines |
1122 | call serialize_pline(mark%pline,bufline,length) |
1123 | @@ -3209,16 +2924,17 @@ MODULE fdf |
1124 | call convert_string_to_array_of_chars(bufline,buffer(init:final)) |
1125 | mark => mark%next |
1126 | enddo |
1127 | - end subroutine serialize_fdf_struct |
1128 | + end subroutine fdf_serialize_struct |
1129 | |
1130 | - subroutine recreate_fdf_struct(nlines,bufferFDF) |
1131 | - character, pointer :: bufferFDF(:) |
1132 | - integer, intent(in) :: nlines |
1133 | + subroutine fdf_recreate_struct(bufferFDF) |
1134 | + character(len=1), intent(in) :: bufferFDF(:) |
1135 | |
1136 | character(len=SERIALIZED_LENGTH) bufline |
1137 | type(parsed_line), pointer :: pline |
1138 | - integer :: i, init, final |
1139 | + integer(ip) :: nlines, i, init, final |
1140 | |
1141 | + nlines = size(bufferFDF) / SERIALIZED_LENGTH |
1142 | + |
1143 | do i= 1, nlines |
1144 | init = (i-1)*SERIALIZED_LENGTH+1 |
1145 | final = (i)*SERIALIZED_LENGTH |
1146 | @@ -3227,121 +2943,14 @@ MODULE fdf |
1147 | call recreate_pline(pline,bufline) |
1148 | call fdf_addtoken(pline%line,pline) |
1149 | enddo |
1150 | - ! print *, "Processed: ", file_in%nlines, " lines." |
1151 | |
1152 | - end subroutine recreate_fdf_struct |
1153 | + end subroutine fdf_recreate_struct |
1154 | |
1155 | -!================================================================= |
1156 | -!--------- Obsolete routines ---------------- |
1157 | -! |
1158 | -! The owner node of the input file sends data file to the |
1159 | -! other processes in the MPI communicator. |
1160 | -! |
1161 | -#ifdef HAVE_CLUSTER_IO |
1162 | - SUBROUTINE fdf_sendInput() |
1163 | - use mpi_siesta |
1164 | - implicit none |
1165 | -!--------------------------------------------------------------- Local Variables |
1166 | - character, pointer :: bufferFDF(:) => null() |
1167 | - integer(ip) :: i, j, k, ierr |
1168 | - type(line_dlist), pointer :: mark |
1169 | + ! To enable client-side setting |
1170 | + SUBROUTINE fdf_set_started(status) |
1171 | + logical, intent(in) :: status |
1172 | |
1173 | -!------------------------------------------------------------------------- BEGIN |
1174 | - ALLOCATE(bufferFDF(file_in%nlines*MAX_LENGTH), stat=ierr) |
1175 | - if (ierr .ne. 0) then |
1176 | - call die('FDF module: fdf_sendInput', 'Error allocating bufferFDF', & |
1177 | - THIS_FILE, __LINE__, fdf_err, rc=ierr) |
1178 | - endif |
1179 | - |
1180 | - mark => file_in%first |
1181 | - do i= 1, file_in%nlines*MAX_LENGTH, MAX_LENGTH |
1182 | - bufferFDF(i:i+MAX_LENGTH-1) = s2arr(mark%str) |
1183 | - mark => mark%next |
1184 | - enddo |
1185 | - |
1186 | - call MPI_Bcast(file_in%nlines, 1, & |
1187 | - MPI_INTEGER, rank, MPI_COMM_WORLD, ierr) |
1188 | - |
1189 | - if (ierr .ne. MPI_SUCCESS) then |
1190 | - call die('FDF module: fdf_sendInput', 'Error Broadcasting nlines.' // & |
1191 | - 'Terminating.', THIS_FILE, __LINE__, fdf_err, rc=ierr) |
1192 | - endif |
1193 | - |
1194 | - call MPI_Bcast(bufferFDF, file_in%nlines*MAX_LENGTH, & |
1195 | - MPI_CHARACTER, rank, MPI_COMM_WORLD, ierr) |
1196 | - if (ierr .ne. MPI_SUCCESS) then |
1197 | - call die('FDF module: fdf_sendInput', 'Error Broadcasting bufferFDF.' // & |
1198 | - 'Terminating.', THIS_FILE, __LINE__, fdf_err, rc=ierr) |
1199 | - endif |
1200 | - |
1201 | - DEALLOCATE(bufferFDF) |
1202 | - |
1203 | - RETURN |
1204 | -!--------------------------------------------------------------------------- END |
1205 | - END SUBROUTINE fdf_sendInput |
1206 | -#endif |
1207 | - |
1208 | -! |
1209 | -! Remaining nodes recieve the input data file from the owner |
1210 | -! of the FDF file inside the cluster. |
1211 | -! |
1212 | -#ifdef HAVE_CLUSTER_IO |
1213 | - SUBROUTINE fdf_recvInput(root, filein, fileinTmp) |
1214 | - use mpi_siesta |
1215 | - implicit none |
1216 | -!--------------------------------------------------------------- Input Variables |
1217 | - character(*) :: filein, fileinTmp |
1218 | - integer(ip) :: root |
1219 | + fdf_started = status |
1220 | + end SUBROUTINE fdf_set_started |
1221 | |
1222 | -!--------------------------------------------------------------- Local Variables |
1223 | - integer(ip) :: i, j, lun, ierr, nlines |
1224 | - character, pointer :: bufferFDF(:) => null() |
1225 | - character(len=10) :: fmt |
1226 | -!------------------------------------------------------------------------- BEGIN |
1227 | - call MPI_Bcast(nlines, 1, & |
1228 | - MPI_INTEGER, root, MPI_COMM_WORLD, ierr) |
1229 | - if (ierr .ne. MPI_SUCCESS) then |
1230 | - call die('FDF module: fdf_recvInput', 'Error Broadcasting nlines.' // & |
1231 | - 'Terminating.', THIS_FILE, __LINE__, fdf_err, rc=ierr) |
1232 | - endif |
1233 | - |
1234 | - ALLOCATE(bufferFDF(nlines*MAX_LENGTH), stat=ierr) |
1235 | - if (ierr .ne. 0) then |
1236 | - call die('FDF module: fdf_recvInput', 'Error allocating bufferFDF', & |
1237 | - THIS_FILE, __LINE__, fdf_err, rc=ierr) |
1238 | - endif |
1239 | - |
1240 | - call MPI_Bcast(bufferFDF, nlines*MAX_LENGTH, & |
1241 | - MPI_CHARACTER, root, MPI_COMM_WORLD, ierr) |
1242 | - if (ierr .ne. MPI_SUCCESS) then |
1243 | - call die('FDF module: fdf_recvInput', 'Error Broadcasting bufferFDF.' // & |
1244 | - 'Terminating.', THIS_FILE, __LINE__, fdf_err, rc=ierr) |
1245 | - endif |
1246 | - |
1247 | - call io_assign(lun) |
1248 | - fileinTmp = TRIM(filein) // '.' // i2s(rank) |
1249 | - open(unit=lun, file=fileinTmp, form='formatted', & |
1250 | - status='unknown') |
1251 | - |
1252 | - if (MAX_LENGTH.lt.10) then |
1253 | - write(fmt,"(a1,I1,a2)") "(", MAX_LENGTH, "a)" |
1254 | - else if (MAX_LENGTH.lt.100) then |
1255 | - write(fmt,"(a1,I2,a2)") "(", MAX_LENGTH, "a)" |
1256 | - else if (MAX_LENGTH.lt.1000) then |
1257 | - write(fmt,"(a1,I3,a2)") "(", MAX_LENGTH, "a)" |
1258 | - else |
1259 | - write(fmt,"(a1,I4,a2)") "(", MAX_LENGTH, "a)" |
1260 | - endif |
1261 | - |
1262 | - do i= 1, nlines*MAX_LENGTH, MAX_LENGTH |
1263 | - write(lun,fmt) (bufferFDF(j),j=i,i+MAX_LENGTH-1) |
1264 | - enddo |
1265 | - call io_close(lun) |
1266 | - |
1267 | - DEALLOCATE(bufferFDF) |
1268 | - |
1269 | - RETURN |
1270 | -!--------------------------------------------------------------------------- END |
1271 | - END SUBROUTINE fdf_recvInput |
1272 | -#endif |
1273 | END MODULE fdf |
1274 | diff --git a/src/utils.F90 b/src/utils.F90 |
1275 | index ee66e11..110ffa2 100644 |
1276 | --- a/src/utils.F90 |
1277 | +++ b/src/utils.F90 |
1278 | @@ -408,12 +408,9 @@ MODULE utils |
1279 | |
1280 | ! |
1281 | ! Die routine (Abort/Terminate program) |
1282 | -! MPI-awareness (MPI_Abort) |
1283 | ! |
1284 | SUBROUTINE die(routine, msg, file, line, unit, rc, cline) |
1285 | -#if defined(HAVE_CLUSTER_IO) || defined(HAVE_BLOCKING_IO) |
1286 | - use mpi_siesta |
1287 | -#endif |
1288 | + |
1289 | implicit none |
1290 | !--------------------------------------------------------------- Input Variables |
1291 | character(len=*), intent(in) :: routine, msg |
1292 | @@ -422,10 +419,6 @@ MODULE utils |
1293 | |
1294 | !--------------------------------------------------------------- Local Variables |
1295 | integer(ip) :: die_unit |
1296 | -#if defined(HAVE_CLUSTER_IO) || defined(HAVE_BLOCKING_IO) |
1297 | - integer(ip) :: ierr, rc2 |
1298 | -#endif |
1299 | - |
1300 | !------------------------------------------------------------------------- BEGIN |
1301 | if (PRESENT(unit)) then |
1302 | die_unit = unit |
1303 | @@ -447,17 +440,8 @@ MODULE utils |
1304 | write(die_unit,'(a)') 'Stopping Program' |
1305 | endif |
1306 | |
1307 | -#if defined(HAVE_CLUSTER_IO) || defined(HAVE_BLOCKING_IO) |
1308 | - if (.not. PRESENT(rc)) then |
1309 | - rc2 = 0 |
1310 | - else |
1311 | - rc2 = rc |
1312 | - endif |
1313 | - |
1314 | - call MPI_Abort(MPI_COMM_WORLD, rc2, ierr) |
1315 | -#else |
1316 | + ! Replace this by a call to a proper handler |
1317 | STOP 'Stopping Program' |
1318 | -#endif |
1319 | !--------------------------------------------------------------------------- END |
1320 | END SUBROUTINE die |
1321 |
This is a welcome change, since it will make the use of LibFDF easier and simplify its maintenance.
If I understand well, I should remove MPI support from the build system once you're done with the merge, right?