Statistiques
| Branche: | Tag: | Révision :

dockonsurf / modules / dos_input.py @ 9cd032cf

Historique | Voir | Annoter | Télécharger (34 ko)

1
"""Functions to deal with DockOnSurf input files.
2

3
Functions
4
try_command:Tries to run a command and logs its exceptions (expected and not).
5
str2lst: Converts a string of integers, and groups of them, to a list of lists.
6
check_expect_val: Checks whether the value of an option has an adequate value.
7
read_input: Sets up the calculation by reading the parameters from input file.
8
get_run_type: Gets 'run_type' value and checks that its value is acceptable.
9
get_code: Gets 'code' value and checks that its value is acceptable.
10
get_batch_q_sys: Gets 'batch_q_sys' value and checks that its value is
11
acceptable.
12
get_relaunch_err: Gets 'relaunch_err' value and checks that its value is
13
acceptable.
14
get_max_jobs: Gets 'max_jobs' value and checks that its value is acceptable.
15
get_special_atoms: Gets 'special_atoms' value and checks that its value is
16
acceptable.
17
get_isol_inp_file: Gets 'isol_inp_file' value and checks that its value is
18
acceptable.
19
get_cluster_magns: Gets 'cluster_magns' value and checks that its value is
20
acceptable.
21
get_num_conformers: Gets 'num_conformers' value and checks that its value is
22
acceptable.
23
get_num_prom_cand: Gets 'num_prom_cand' value and checks that its value is
24
acceptable.
25
get_iso_rmsd: Gets 'iso_rmsd' value and checks that its value is acceptable.
26
get_pre_opt: Gets 'pre_opt' value and checks that its value is acceptable.
27
get_screen_inp_file: Gets 'screen_inp_file' value and checks that its value is
28
acceptable.
29
get_sites: Gets 'sites' value and checks that its value is acceptable.
30
get_molec_ctrs: Gets 'molec_ctrs' value and checks that its value is
31
acceptable.
32
get_try_disso: Gets 'try_disso' value and checks that its value is acceptable.
33
get_pts_per_angle: Gets 'pts_per_angle' value and checks that its value is
34
acceptable.
35
get_coll_thrsld: Gets 'coll_thrsld' value and checks that its value is
36
acceptable.
37
get_screen_rmsd: Gets 'screen_rmsd' value and checks that its value is
38
acceptable.
39
get_coll_bottom_z: Gets 'coll_bottom_z' value and checks that its value is
40
acceptable.
41
get_refine_inp_file: Gets 'refine_inp_file' value and checks that its value is
42
acceptable.
43
get_energy_cutoff: Gets 'energy_cutoff' value and checks that its value is
44
acceptable.
45
"""
46
import os.path
47
import logging
48
from configparser import ConfigParser, NoSectionError, NoOptionError, \
49
    MissingSectionHeaderError, DuplicateOptionError
50
import numpy as np
51
from modules.utilities import try_command
52

    
53
logger = logging.getLogger('DockOnSurf')
54

    
55
dos_inp = ConfigParser(inline_comment_prefixes='#',
56
                       empty_lines_in_values=False)
57

    
58
new_answers = {'n': False, 'none': False, 'nay': False,
59
               'y': True, '': True, 'aye': True, 'sure': True}
60
for answer, val in new_answers.items():
61
    dos_inp.BOOLEAN_STATES[answer] = val  # TODO Check value 0
62
turn_false_answers = [answer for answer in dos_inp.BOOLEAN_STATES
63
                      if dos_inp.BOOLEAN_STATES[answer] is False]
64
turn_true_answers = [answer for answer in dos_inp.BOOLEAN_STATES
65
                     if dos_inp.BOOLEAN_STATES[answer]]
66

    
67
no_sect_err = "Section '%s' not found on input file"
68
no_opt_err = "Option '%s' not found on section '%s'"
69
num_error = "'%s' value must be a %s"
70

    
71

    
72
def str2lst(cmplx_str, func=int):  # TODO: enable deeper level of nested lists
73
    # TODO Treat all-enclosing parenthesis as a list instead of list of lists.
74
    """Converts a string of integers, and groups of them, to a list.
75

76
    Keyword arguments:
77
    @param cmplx_str: str, string of integers (or floats) and groups of them
78
    enclosed by parentheses-like characters.
79
    - Group enclosers: '()' '[]' and '{}'.
80
    - Separators: ',' ';' and ' '.
81
    - Nested groups are not allowed: '3 ((6 7) 8) 4'.
82
    @param func: either to use int or float
83

84
    @return list, list of integers (or floats), or list of integers (or floats)
85
    in the case they were grouped. First, the singlets are placed, and then the
86
    groups in input order.
87

88
    eg. '128,(135 138;141] 87 {45, 68}' -> [128, 87, [135, 138, 141], [45, 68]]
89
    """
90

    
91
    # Checks
92
    error_msg = "Function argument should be a str, sequence of " \
93
                "numbers separated by ',' ';' or ' '." \
94
                "\nThey can be grouped in parentheses-like enclosers: '()', " \
95
                "'[]' or {}. Nested groups are not allowed. \n" \
96
                "eg. 128,(135 138;141) 87 {45, 68}"
97
    cmplx_str = try_command(cmplx_str.replace, [(AttributeError, error_msg)],
98
                            ',', ' ')
99

    
100
    cmplx_str = cmplx_str.replace(';', ' ').replace('[', '(').replace(
101
        ']', ')').replace('{', '(').replace('}', ')')
102

    
103
    try_command(list, [(ValueError, error_msg)], map(func, cmplx_str.replace(
104
        ')', '').replace('(', '').split()))
105

    
106
    deepness = 0
107
    for el in cmplx_str.split():
108
        if '(' in el:
109
            deepness += 1
110
        if ')' in el:
111
            deepness += -1
112
        if deepness > 1 or deepness < 0:
113
            logger.error(error_msg)
114
            raise ValueError(error_msg)
115

    
116
    init_list = cmplx_str.split()
117
    start_group = []
118
    end_group = []
119
    for i, element in enumerate(init_list):
120
        if '(' in element:
121
            start_group.append(i)
122
            init_list[i] = element.replace('(', '')
123
        if ')' in element:
124
            end_group.append(i)
125
            init_list[i] = element.replace(')', '')
126

    
127
    init_list = list(map(func, init_list))
128

    
129
    new_list = []
130
    for start_el, end_el in zip(start_group, end_group):
131
        new_list.append(init_list[start_el:end_el + 1])
132

    
133
    for v in new_list:
134
        for el in v:
135
            init_list.remove(el)
136
    return init_list + new_list
137

    
138

    
139
def check_expect_val(value, expect_vals, err_msg=None):
140
    """Checks whether an option lies within its expected values.
141

142
    Keyword arguments:
143
    @param value: The variable to check if its value lies within the expected
144
    ones
145
    @param expect_vals: list, list of values allowed for the present option.
146
    @param err_msg: The error message to be prompted in both log and screen.
147
    @raise ValueError: if the value is not among the expected ones.
148
    @return True if the value is among the expected ones.
149
    """
150
    if err_msg is None:
151
        err_msg = f"'{value}' is not an adequate value.\n" \
152
                  f"Adequate values: {expect_vals}"
153
    if not any([exp_val == value for exp_val in expect_vals]):
154
        logger.error(err_msg)
155
        raise ValueError(err_msg)
156

    
157
    return True
158

    
159

    
160
def check_inp_file(inp_file, code):
161
    if code == 'cp2k':
162
        from pycp2k import CP2K
163
        cp2k = CP2K()
164
        try_command(cp2k.parse,
165
                    [(UnboundLocalError, "Invalid CP2K input file")], inp_file)
166

    
167

    
168
# Global
169

    
170
def get_run_type():
171
    isolated, screening, refinement = (False, False, False)
172
    run_type_vals = ['isolated', 'screening', 'refinement', 'adsorption',
173
                     'full']
174
    run_types = dos_inp.get('Global', 'run_type').split()
175
    for run_type in run_types:
176
        check_expect_val(run_type.lower(), run_type_vals)
177
        if 'isol' in run_type.lower():
178
            isolated = True
179
        if 'screen' in run_type.lower():
180
            screening = True
181
        if 'refine' in run_type.lower():
182
            refinement = True
183
        if 'adsor' in run_type.lower():
184
            screening, refinement = (True, True)
185
        if 'full' in run_type.lower():
186
            isolated, screening, refinement = (True, True, True)
187

    
188
    return isolated, screening, refinement
189

    
190

    
191
def get_code():
192
    code_vals = ['cp2k']
193
    check_expect_val(dos_inp.get('Global', 'code').lower(), code_vals)
194
    code = dos_inp.get('Global', 'code').lower()
195
    return code
196

    
197

    
198
def get_batch_q_sys():
199
    batch_q_sys_vals = ['sge', 'lsf', 'irene', 'local'] + turn_false_answers
200
    check_expect_val(dos_inp.get('Global', 'batch_q_sys').lower(),
201
                     batch_q_sys_vals)
202
    batch_q_sys = dos_inp.get('Global', 'batch_q_sys').lower()
203
    if batch_q_sys.lower() in turn_false_answers:
204
        return False
205
    else:
206
        return batch_q_sys
207

    
208

    
209
def get_subm_script():
210
    subm_script = dos_inp.get('Global', 'subm_script', fallback=False)
211
    if subm_script and not os.path.isfile(subm_script):
212
        logger.error(f'File {subm_script} not found.')
213
        raise FileNotFoundError(f'File {subm_script} not found')
214
    return subm_script
215

    
216

    
217
def get_project_name():
218
    project_name = dos_inp.get('Global', 'project_name', fallback='')
219
    return project_name
220

    
221

    
222
def get_relaunch_err():
223
    relaunch_err_vals = ['geo_not_conv']
224
    relaunch_err = dos_inp.get('Global', 'relaunch_err',
225
                               fallback="False")
226
    if relaunch_err.lower() in turn_false_answers:
227
        return False
228
    else:
229
        check_expect_val(relaunch_err.lower(), relaunch_err_vals)
230
    return relaunch_err
231

    
232

    
233
def get_max_jobs():
234
    import re
235
    err_msg = "'max_jobs' must be a list of, number plus 'p', 'q' or 'r', or " \
236
              "a combination of them without repeating letters.\n" \
237
              "eg: '2r 3p 4pr', '5q' or '3r 3p'"
238
    max_jobs_str = dos_inp.get('Global', 'max_jobs', fallback="inf").lower()
239
    str_vals = ["r", "p", "q", "rp", "rq", "pr", "qr"]
240
    max_jobs = {"r": np.inf, "p": np.inf, "rp": np.inf}
241
    if "inf" == max_jobs_str:
242
        return {"r": np.inf, "p": np.inf, "rp": np.inf}
243
    # Iterate over the number of requirements:
244
    for req in max_jobs_str.split():
245
        # Split numbers from letters into a list
246
        req_parts = re.findall(r'[a-z]+|\d+', req)
247
        if len(req_parts) != 2:
248
            logger.error(err_msg)
249
            raise ValueError(err_msg)
250
        if req_parts[0].isdecimal():
251
            req_parts[1] = req_parts[1].replace('q', 'p').replace('pr', 'rp')
252
            if req_parts[1] in str_vals and max_jobs[req_parts[1]] == np.inf:
253
                max_jobs[req_parts[1]] = int(req_parts[0])
254
        elif req_parts[1].isdecimal():
255
            req_parts[0] = req_parts[0].replace('q', 'p').replace('pr', 'rp')
256
            if req_parts[0] in str_vals and max_jobs[req_parts[0]] == np.inf:
257
                max_jobs[req_parts[0]] = int(req_parts[1])
258
        else:
259
            logger.error(err_msg)
260
            raise ValueError(err_msg)
261

    
262
    return max_jobs
263

    
264

    
265
def get_special_atoms():
266
    from ase.data import chemical_symbols
267

    
268
    spec_at_err = '\'special_atoms\' does not have an adequate format.\n' \
269
                  'Adequate format: (Fe1 Fe) (O1 O)'
270
    special_atoms = dos_inp.get('Global', 'special_atoms', fallback="False")
271
    if special_atoms.lower() in turn_false_answers:
272
        special_atoms = []
273
    else:
274
        # Converts the string into a list of tuples
275
        lst_tple = [tuple(pair.replace("(", "").split()) for pair in
276
                    special_atoms.split(")")[:-1]]
277
        if len(lst_tple) == 0:
278
            logger.error(spec_at_err)
279
            raise ValueError(spec_at_err)
280
        for i, tup in enumerate(lst_tple):
281
            if not isinstance(tup, tuple) or len(tup) != 2:
282
                logger.error(spec_at_err)
283
                raise ValueError(spec_at_err)
284
            if tup[1].capitalize() not in chemical_symbols:
285
                elem_err = "The second element of the couple should be an " \
286
                           "actual element of the periodic table"
287
                logger.error(elem_err)
288
                raise ValueError(elem_err)
289
            if tup[0].capitalize() in chemical_symbols:
290
                elem_err = "The first element of the couple is already an " \
291
                           "actual element of the periodic table, "
292
                logger.error(elem_err)
293
                raise ValueError(elem_err)
294
            for j, tup2 in enumerate(lst_tple):
295
                if j <= i:
296
                    continue
297
                if tup2[0] == tup[0]:
298
                    label_err = f'You have specified the label {tup[0]} to ' \
299
                                f'more than one special atom'
300
                    logger.error(label_err)
301
                    raise ValueError(label_err)
302
        special_atoms = lst_tple
303
    return special_atoms
304

    
305

    
306
# Isolated
307

    
308
def get_isol_inp_file():
309
    isol_inp_file = dos_inp.get('Isolated', 'isol_inp_file')
310
    if not os.path.isfile(isol_inp_file):
311
        logger.error(f'File {isol_inp_file} not found.')
312
        raise FileNotFoundError(f'File {isol_inp_file} not found')
313
    return isol_inp_file
314

    
315

    
316
def get_molec_file():
317
    molec_file = dos_inp.get('Isolated', 'molec_file')
318
    if not os.path.isfile(molec_file):
319
        logger.error(f'File {molec_file} not found.')
320
        raise FileNotFoundError(f'File {molec_file} not found')
321
    return molec_file
322

    
323

    
324
def get_num_conformers():
325
    err_msg = num_error % ('num_conformers', 'positive integer')
326
    num_conformers = try_command(dos_inp.getint, [(ValueError, err_msg)],
327
                                 'Isolated', 'num_conformers', fallback=100)
328
    if num_conformers < 1:
329
        logger.error(err_msg)
330
        raise ValueError(err_msg)
331
    return num_conformers
332

    
333

    
334
def get_pre_opt():
335
    pre_opt_vals = ['uff', 'mmff'] + turn_false_answers
336
    check_expect_val(dos_inp.get('Isolated', 'pre_opt').lower(), pre_opt_vals)
337
    pre_opt = dos_inp.get('Isolated', 'pre_opt').lower()
338
    if pre_opt in turn_false_answers:
339
        return False
340
    else:
341
        return pre_opt
342

    
343

    
344
# Screening
345

    
346
def get_screen_inp_file():
347
    screen_inp_file = dos_inp.get('Screening', 'screen_inp_file')
348
    if not os.path.isfile(screen_inp_file):
349
        logger.error(f'File {screen_inp_file} not found.')
350
        raise FileNotFoundError(f'File {screen_inp_file} not found')
351
    return screen_inp_file
352

    
353

    
354
def get_surf_file():
355
    surf_file = dos_inp.get('Screening', 'surf_file')
356
    if not os.path.isfile(surf_file):
357
        logger.error(f'File {surf_file} not found.')
358
        raise FileNotFoundError(f'File {surf_file} not found')
359
    return surf_file
360

    
361

    
362
def get_sites():
363
    err_msg = 'The value of sites should be a list of atom numbers ' \
364
              '(ie. positive integers) or groups of atom numbers ' \
365
              'grouped by parentheses-like enclosers. \n' \
366
              'eg. 128,(135 138;141) 87 {45, 68}'
367
    # Convert the string into a list of lists
368
    sites = try_command(str2lst,
369
                        [(ValueError, err_msg), (AttributeError, err_msg)],
370
                        dos_inp.get('Screening', 'sites'))
371
    # Check all elements of the list (of lists) are positive integers
372
    for site in sites:
373
        if type(site) is list:
374
            for atom in site:
375
                if atom < 0:
376
                    logger.error(err_msg)
377
                    raise ValueError(err_msg)
378
        elif type(site) is int:
379
            if site < 0:
380
                logger.error(err_msg)
381
                raise ValueError(err_msg)
382
        else:
383
            logger.error(err_msg)
384
            raise ValueError(err_msg)
385

    
386
    return sites
387

    
388

    
389
def get_surf_ctrs2():
390
    err_msg = 'The value of surf_ctrs2 should be a list of atom numbers ' \
391
              '(ie. positive integers) or groups of atom numbers ' \
392
              'grouped by parentheses-like enclosers. \n' \
393
              'eg. 128,(135 138;141) 87 {45, 68}'
394
    # Convert the string into a list of lists
395
    surf_ctrs2 = try_command(str2lst,
396
                             [(ValueError, err_msg), (AttributeError, err_msg)],
397
                             dos_inp.get('Screening', 'surf_ctrs2'))
398
    # Check all elements of the list (of lists) are positive integers
399
    for ctr in surf_ctrs2:
400
        if type(ctr) is list:
401
            for atom in ctr:
402
                if atom < 0:
403
                    logger.error(err_msg)
404
                    raise ValueError(err_msg)
405
        elif type(ctr) is int:
406
            if ctr < 0:
407
                logger.error(err_msg)
408
                raise ValueError(err_msg)
409
        else:
410
            logger.error(err_msg)
411
            raise ValueError(err_msg)
412

    
413
    return surf_ctrs2
414

    
415

    
416
def get_molec_ctrs():
417
    err_msg = 'The value of molec_ctrs should be a list of atom' \
418
              ' numbers (ie. positive integers) or groups of atom ' \
419
              'numbers enclosed by parentheses-like characters. \n' \
420
              'eg. 128,(135 138;141) 87 {45, 68}'
421
    # Convert the string into a list of lists
422
    molec_ctrs = try_command(str2lst,
423
                             [(ValueError, err_msg),
424
                              (AttributeError, err_msg)],
425
                             dos_inp.get('Screening', 'molec_ctrs'))
426
    # Check all elements of the list (of lists) are positive integers
427
    for ctr in molec_ctrs:
428
        if isinstance(ctr, list):
429
            for atom in ctr:
430
                if atom < 0:
431
                    logger.error(err_msg)
432
                    raise ValueError(err_msg)
433
        elif isinstance(ctr, int):
434
            if ctr < 0:
435
                logger.error(err_msg)
436
                raise ValueError(err_msg)
437
        else:
438
            logger.error(err_msg)
439
            raise ValueError(err_msg)
440

    
441
    return molec_ctrs
442

    
443

    
444
def get_molec_ctrs2():
445
    err_msg = 'The value of molec_ctrs2 should be a list of atom ' \
446
              'numbers (ie. positive integers) or groups of atom ' \
447
              'numbers enclosed by parentheses-like characters. \n' \
448
              'eg. 128,(135 138;141) 87 {45, 68}'
449
    # Convert the string into a list of lists
450
    molec_ctrs2 = try_command(str2lst, [(ValueError, err_msg),
451
                                        (AttributeError, err_msg)],
452
                              dos_inp.get('Screening', 'molec_ctrs2'))
453

    
454
    # Check all elements of the list (of lists) are positive integers
455
    for ctr in molec_ctrs2:
456
        if isinstance(ctr, list):
457
            for atom in ctr:
458
                if atom < 0:
459
                    logger.error(err_msg)
460
                    raise ValueError(err_msg)
461
        elif isinstance(ctr, int):
462
            if ctr < 0:
463
                logger.error(err_msg)
464
                raise ValueError(err_msg)
465
        else:
466
            logger.error(err_msg)
467
            raise ValueError(err_msg)
468

    
469
    return molec_ctrs2
470

    
471

    
472
def get_molec_ctrs3():
473
    err_msg = 'The value of molec_ctrs3 should be a list of atom ' \
474
              'numbers (ie. positive integers) or groups of atom ' \
475
              'numbers enclosed by parentheses-like characters. \n' \
476
              'eg. 128,(135 138;141) 87 {45, 68}'
477
    # Convert the string into a list of lists
478
    molec_ctrs3 = try_command(str2lst, [(ValueError, err_msg),
479
                                        (AttributeError, err_msg)],
480
                              dos_inp.get('Screening', 'molec_ctrs3'))
481

    
482
    # Check all elements of the list (of lists) are positive integers
483
    for ctr in molec_ctrs3:
484
        if isinstance(ctr, list):
485
            for atom in ctr:
486
                if atom < 0:
487
                    logger.error(err_msg)
488
                    raise ValueError(err_msg)
489
        elif isinstance(ctr, int):
490
            if ctr < 0:
491
                logger.error(err_msg)
492
                raise ValueError(err_msg)
493
        else:
494
            logger.error(err_msg)
495
            raise ValueError(err_msg)
496

    
497
    return molec_ctrs3
498

    
499

    
500
def get_max_helic_angle():
501
    err_msg = "'max_helic_angle' must be a positive number in degrees"
502
    max_helic_angle = try_command(dos_inp.getfloat, [(ValueError, err_msg)],
503
                                  'Screening', 'max_helic_angle',
504
                                  fallback=180.0)
505
    if max_helic_angle < 0:
506
        logger.error(err_msg)
507
        raise ValueError(err_msg)
508

    
509
    return max_helic_angle
510

    
511

    
512
def get_select_magns():
513
    select_magns_vals = ['energy', 'moi']
514
    select_magns_str = dos_inp.get('Screening', 'select_magns',
515
                                   fallback='energy')
516
    select_magns_str.replace(',', ' ').replace(';', ' ')
517
    select_magns = select_magns_str.split(' ')
518
    select_magns = [m.lower() for m in select_magns]
519
    for m in select_magns:
520
        check_expect_val(m, select_magns_vals)
521
    return select_magns
522

    
523

    
524
def get_confs_per_magn():
525
    err_msg = num_error % ('confs_per_magn', 'positive integer')
526
    confs_per_magn = try_command(dos_inp.getint, [(ValueError, err_msg)],
527
                                 'Screening', 'confs_per_magn', fallback=2)
528
    if confs_per_magn <= 0:
529
        logger.error(err_msg)
530
        raise ValueError(err_msg)
531
    return confs_per_magn
532

    
533

    
534
def get_pbc_cell():
535
    err = "'pbc_cell' must be either 3 vectors of size 3 or False"
536
    pbc_cell_str = dos_inp.get('Screening', 'pbc_cell', fallback="False")
537
    if pbc_cell_str.lower() in turn_false_answers:
538
        return False
539
    else:
540
        pbc_cell = np.array(try_command(str2lst, [(ValueError, err)],
541
                                        pbc_cell_str, float))
542
        if pbc_cell.shape != (3, 3):
543
            logger.error(err)
544
            raise ValueError(err)
545
        return pbc_cell
546

    
547

    
548
def get_surf_norm_vect():
549
    err = "'surf_norm_vect' must be a 3 component vector, 'x', 'y' or 'z', " \
550
          "'auto' or 'asann'."
551
    cart_axes = {'x': [1.0, 0.0, 0.0], '-x': [-1.0, 0.0, 0.0],
552
                 'y': [0.0, 1.0, 0.0], '-y': [0.0, -1.0, 0.0],
553
                 'z': [0.0, 0.0, 1.0], '-z': [0.0, 0.0, -1.0]}
554
    surf_norm_vect_str = dos_inp.get('Screening', 'surf_norm_vect',
555
                                     fallback="auto").lower()
556
    if surf_norm_vect_str == "asann" or surf_norm_vect_str == "auto":
557
        return 'auto'
558
    if surf_norm_vect_str in cart_axes:
559
        return np.array(cart_axes[surf_norm_vect_str])
560
    surf_norm_vect = try_command(str2lst, [(ValueError, err)],
561
                                 surf_norm_vect_str, float)
562
    if len(surf_norm_vect) != 3:
563
        logger.error(err)
564
        raise ValueError(err)
565

    
566
    return np.array(surf_norm_vect) / np.linalg.norm(surf_norm_vect)
567

    
568

    
569
def get_adsorption_height():
570
    err_msg = num_error % ('adsorption_height', 'positive number')
571
    ads_height = try_command(dos_inp.getfloat, [(ValueError, err_msg)],
572
                             'Screening', 'adsorption_height', fallback=2.5)
573
    if ads_height <= 0:
574
        logger.error(err_msg)
575
        raise ValueError(err_msg)
576
    return ads_height
577

    
578

    
579
def get_set_angles():
580
    set_vals = ['euler', 'internal']
581
    check_expect_val(dos_inp.get('Screening', 'set_angles').lower(), set_vals)
582
    set_angles = dos_inp.get('Screening', 'set_angles',
583
                             fallback='euler').lower()
584
    return set_angles
585

    
586

    
587
def get_pts_per_angle():
588
    err_msg = num_error % ('sample_points_per_angle', 'positive integer')
589
    pts_per_angle = try_command(dos_inp.getint,
590
                                [(ValueError, err_msg)],
591
                                'Screening', 'sample_points_per_angle',
592
                                fallback=3)
593
    if pts_per_angle <= 0:
594
        logger.error(err_msg)
595
        raise ValueError(err_msg)
596
    return pts_per_angle
597

    
598

    
599
def get_max_structures():
600
    err_msg = num_error % ('max_structures', 'positive integer')
601
    max_structs = dos_inp.get('Screening', 'max_structures', fallback="False")
602
    if max_structs.lower() in turn_false_answers:
603
        return np.inf
604
    if try_command(int, [(ValueError, err_msg)], max_structs) <= 0:
605
        logger.error(err_msg)
606
        raise ValueError(err_msg)
607
    return int(max_structs)
608

    
609

    
610
def get_coll_thrsld():
611
    err_msg = num_error % ('collision_threshold', 'positive number')
612
    coll_thrsld_str = dos_inp.get('Screening', 'collision_threshold',
613
                                  fallback="1.3")
614
    if coll_thrsld_str.lower() in turn_false_answers:
615
        return False
616
    coll_thrsld = try_command(float, [(ValueError, err_msg)], coll_thrsld_str)
617

    
618
    if coll_thrsld <= 0:
619
        logger.error(err_msg)
620
        raise ValueError(err_msg)
621

    
622
    return coll_thrsld
623

    
624

    
625
def get_min_coll_height(norm_vect):
626
    err_msg = num_error % ('min_coll_height', 'decimal number')
627
    min_coll_height = dos_inp.get('Screening', 'min_coll_height',
628
                                  fallback="False")
629
    if min_coll_height.lower() in turn_false_answers:
630
        return False
631
    min_coll_height = try_command(float, [(ValueError, err_msg)],
632
                                  min_coll_height)
633
    cart_axes = [[1.0, 0.0, 0.0], [0.0, 1.0, 0.0], [0.0, 0.0, 1.0],
634
                 [-1.0, 0.0, 0.0], [0.0, -1.0, 0.0], [0.0, 0.0, -1.0]]
635
    err_msg = "'min_coll_height' option is only implemented for " \
636
              "'surf_norm_vect' to be one of the x, y or z axes. "
637
    if not isinstance(norm_vect, str) or norm_vect != 'auto':
638
        check_expect_val(norm_vect.tolist(), cart_axes, err_msg)
639
    return min_coll_height
640

    
641

    
642
def get_exclude_ads_ctr():
643
    err_msg = "exclude_ads_ctr must have a boolean value."
644
    exclude_ads_ctr = try_command(dos_inp.getboolean, [(ValueError, err_msg)],
645
                                  "Screening", "exclude_ads_ctr",
646
                                  fallback=False)
647
    return exclude_ads_ctr
648

    
649

    
650
def get_H_donor(spec_atoms):
651
    from ase.data import chemical_symbols
652
    err_msg = "The value of 'h_donor' must be either False, a chemical symbol " \
653
              "or an atom index"
654
    h_donor_str = dos_inp.get('Screening', 'h_donor', fallback="False")
655
    h_donor = []
656
    if h_donor_str.lower() in turn_false_answers:
657
        return False
658
    err = False
659
    for el in h_donor_str.split():
660
        try:
661
            h_donor.append(int(el))
662
        except ValueError:
663
            if el not in chemical_symbols + [nw_sym for pairs in spec_atoms
664
                                             for nw_sym in pairs]:
665
                err = True
666
            else:
667
                h_donor.append(el)
668
        finally:
669
            if err:
670
                logger.error(err_msg)
671
                ValueError(err_msg)
672
    return h_donor
673

    
674

    
675
def get_H_acceptor(spec_atoms):
676
    from ase.data import chemical_symbols
677
    err_msg = "The value of 'h_acceptor' must be either 'all', a chemical " \
678
              "symbol or an atom index"
679
    h_acceptor_str = dos_inp.get('Screening', 'h_acceptor', fallback="all")
680
    if h_acceptor_str.lower() == "all":
681
        return "all"
682
    h_acceptor = []
683
    err = False
684
    for el in h_acceptor_str.split():
685
        try:
686
            h_acceptor.append(int(el))
687
        except ValueError:
688
            if el not in chemical_symbols + [nw_sym for pairs in spec_atoms
689
                                             for nw_sym in pairs]:
690
                err = True
691
            else:
692
                h_acceptor.append(el)
693
        finally:
694
            if err:
695
                logger.error(err_msg)
696
                ValueError(err_msg)
697
    return h_acceptor
698

    
699

    
700
def get_use_molec_file():
701
    use_molec_file = dos_inp.get('Screening', 'use_molec_file',
702
                                 fallback='False')
703
    if use_molec_file.lower() in turn_false_answers:
704
        return False
705
    if not os.path.isfile(use_molec_file):
706
        logger.error(f'File {use_molec_file} not found.')
707
        raise FileNotFoundError(f'File {use_molec_file} not found')
708

    
709
    return use_molec_file
710

    
711

    
712
# Refinement
713

    
714
def get_refine_inp_file():
715
    refine_inp_file = dos_inp.get('Refinement', 'refine_inp_file')
716
    if not os.path.isfile(refine_inp_file):
717
        logger.error(f'File {refine_inp_file} not found.')
718
        raise FileNotFoundError(f'File {refine_inp_file} not found')
719

    
720
    return refine_inp_file
721

    
722

    
723
def get_energy_cutoff():
724
    err_msg = num_error % ('energy_cutoff', 'positive decimal number')
725
    energy_cutoff = try_command(dos_inp.getfloat,
726
                                [(ValueError, err_msg)],
727
                                'Refinement', 'energy_cutoff', fallback=0.5)
728
    if energy_cutoff < 0:
729
        logger.error(err_msg)
730
        raise ValueError(err_msg)
731
    return energy_cutoff
732

    
733

    
734
# Read input parameters
735

    
736
def read_input(in_file):
737
    err = False
738
    try:
739
        dos_inp.read(in_file)
740
    except MissingSectionHeaderError as e:
741
        logger.error('There are options in the input file without a Section '
742
                     'header.')
743
        err = e
744
    except DuplicateOptionError as e:
745
        logger.error('There is an option in the input file that has been '
746
                     'specified more than once.')
747
        err = e
748
    except Exception as e:
749
        err = e
750
    else:
751
        err = False
752
    finally:
753
        if isinstance(err, BaseException):
754
            raise err
755

    
756
    inp_vars = {}
757

    
758
    # Global
759
    if not dos_inp.has_section('Global'):
760
        logger.error(no_sect_err % 'Global')
761
        raise NoSectionError('Global')
762

    
763
    # Mandatory options
764
    # Checks whether the mandatory options 'run_type', 'code', etc. are present.
765
    glob_mand_opts = ['run_type', 'batch_q_sys']
766
    for opt in glob_mand_opts:
767
        if not dos_inp.has_option('Global', opt):
768
            logger.error(no_opt_err % (opt, 'Global'))
769
            raise NoOptionError(opt, 'Global')
770

    
771
    # Gets which sections are to be carried out
772
    isolated, screening, refinement = get_run_type()
773
    inp_vars['isolated'] = isolated
774
    inp_vars['screening'] = screening
775
    inp_vars['refinement'] = refinement
776
    inp_vars['batch_q_sys'] = get_batch_q_sys()
777

    
778
    # Dependent options:
779
    inp_vars['code'] = get_code()
780
    if inp_vars['batch_q_sys']:
781
        inp_vars['max_jobs'] = get_max_jobs()
782
        if inp_vars['batch_q_sys'] != 'local':
783
            if not dos_inp.has_option('Global', 'subm_script'):
784
                logger.error(no_opt_err % ('subm_script', 'Global'))
785
                raise NoOptionError('subm_script', 'Global')
786
            inp_vars['subm_script'] = get_subm_script()
787

    
788
    # Facultative options (Default/Fallback value present)
789
    inp_vars['project_name'] = get_project_name()
790
    # inp_vars['relaunch_err'] = get_relaunch_err()
791
    inp_vars['special_atoms'] = get_special_atoms()
792

    
793
    # Isolated
794
    if isolated:
795
        if not dos_inp.has_section('Isolated'):
796
            logger.error(no_sect_err % 'Isolated')
797
            raise NoSectionError('Isolated')
798
        # Mandatory options
799
        # Checks whether the mandatory options are present.
800
        iso_mand_opts = ['isol_inp_file', 'molec_file']
801
        for opt in iso_mand_opts:
802
            if not dos_inp.has_option('Isolated', opt):
803
                logger.error(no_opt_err % (opt, 'Isolated'))
804
                raise NoOptionError(opt, 'Isolated')
805
        inp_vars['isol_inp_file'] = get_isol_inp_file()
806
        if 'code ' in inp_vars:
807
            check_inp_file(inp_vars['isol_inp_file'], inp_vars['code'])
808
        inp_vars['molec_file'] = get_molec_file()
809

    
810
        # Facultative options (Default/Fallback value present)
811
        inp_vars['num_conformers'] = get_num_conformers()
812
        # inp_vars['num_prom_cand'] = get_num_prom_cand()
813
        # inp_vars['iso_rmsd'] = get_iso_rmsd()
814
        inp_vars['pre_opt'] = get_pre_opt()
815

    
816
    # Screening
817
    if screening:
818
        if not dos_inp.has_section('Screening'):
819
            logger.error(no_sect_err % 'Screening')
820
            raise NoSectionError('Screening')
821
        # Mandatory options:
822
        # Checks whether the mandatory options are present.
823
        screen_mand_opts = ['screen_inp_file', 'surf_file', 'sites',
824
                            'molec_ctrs']
825
        for opt in screen_mand_opts:
826
            if not dos_inp.has_option('Screening', opt):
827
                logger.error(no_opt_err % (opt, 'Screening'))
828
                raise NoOptionError(opt, 'Screening')
829
        inp_vars['screen_inp_file'] = get_screen_inp_file()
830
        inp_vars['surf_file'] = get_surf_file()
831
        inp_vars['sites'] = get_sites()
832
        inp_vars['molec_ctrs'] = get_molec_ctrs()
833

    
834
        # Facultative options (Default value present)
835
        inp_vars['select_magns'] = get_select_magns()
836
        inp_vars['confs_per_magn'] = get_confs_per_magn()
837
        inp_vars['surf_norm_vect'] = get_surf_norm_vect()
838
        inp_vars['adsorption_height'] = get_adsorption_height()
839
        inp_vars['set_angles'] = get_set_angles()
840
        inp_vars['sample_points_per_angle'] = get_pts_per_angle()
841
        inp_vars['pbc_cell'] = get_pbc_cell()
842
        inp_vars['collision_threshold'] = get_coll_thrsld()
843
        inp_vars['min_coll_height'] = get_min_coll_height(
844
            inp_vars['surf_norm_vect'])
845
        if inp_vars['min_coll_height'] is False \
846
                and inp_vars['collision_threshold'] is False:
847
            logger.warning("Collisions are deactivated: Overlapping of "
848
                           "adsorbate and surface is possible")
849
        inp_vars['exclude_ads_ctr'] = get_exclude_ads_ctr()
850
        inp_vars['h_donor'] = get_H_donor(inp_vars['special_atoms'])
851
        inp_vars['max_structures'] = get_max_structures()
852
        inp_vars['use_molec_file'] = get_use_molec_file()
853

    
854
        # Options depending on the value of others
855
        if inp_vars['set_angles'] == "internal":
856
            internal_opts = ['molec_ctrs2', 'molec_ctrs3', 'surf_ctrs2',
857
                             'max_helic_angle']
858
            for opt in internal_opts:
859
                if not dos_inp.has_option('Screening', opt):
860
                    logger.error(no_opt_err % (opt, 'Screening'))
861
                    raise NoOptionError(opt, 'Screening')
862
            inp_vars['max_helic_angle'] = get_max_helic_angle()
863
            inp_vars['molec_ctrs2'] = get_molec_ctrs2()
864
            inp_vars['molec_ctrs3'] = get_molec_ctrs3()
865
            inp_vars['surf_ctrs2'] = get_surf_ctrs2()
866
            if len(inp_vars["molec_ctrs2"]) != len(inp_vars['molec_ctrs']) \
867
                    or len(inp_vars["molec_ctrs3"]) != \
868
                    len(inp_vars['molec_ctrs']) \
869
                    or len(inp_vars['surf_ctrs2']) != len(inp_vars['sites']):
870
                err = "'molec_ctrs' 'molec_ctrs2' and 'molec_ctrs2' must have " \
871
                      "the same number of indices "
872
                logger.error(err)
873
                raise ValueError(err)
874

    
875
        if inp_vars['h_donor'] is False:
876
            inp_vars['h_acceptor'] = False
877
        else:
878
            inp_vars['h_acceptor'] = get_H_acceptor(inp_vars['special_atoms'])
879

    
880
    # Refinement
881
    if refinement:
882
        if not dos_inp.has_section('Refinement'):
883
            logger.error(no_sect_err % 'Refinement')
884
            raise NoSectionError('Refinement')
885
        # Mandatory options
886
        # Checks whether the mandatory options are present.
887
        ref_mand_opts = ['refine_inp_file']
888
        for opt in ref_mand_opts:
889
            if not dos_inp.has_option('Refinement', opt):
890
                logger.error(no_opt_err % (opt, 'Refinement'))
891
                raise NoOptionError(opt, 'Refinement')
892
        inp_vars['refine_inp_file'] = get_refine_inp_file()
893

    
894
        # Facultative options (Default value present)
895
        inp_vars['energy_cutoff'] = get_energy_cutoff()
896
        # end energy_cutoff
897

    
898
    return_vars_str = "\n\t".join([str(key) + ": " + str(value)
899
                                   for key, value in inp_vars.items()])
900
    logger.info(
901
        f'Correctly read {in_file} parameters: \n\n\t{return_vars_str}\n')
902

    
903
    return inp_vars
904

    
905

    
906
if __name__ == "__main__":
907
    import sys
908

    
909
    print(read_input(sys.argv[1]))