1!
2!  Dalton, a molecular electronic structure program
3!  Copyright (C) by the authors of Dalton.
4!
5!  This program is free software; you can redistribute it and/or
6!  modify it under the terms of the GNU Lesser General Public
7!  License version 2.1 as published by the Free Software Foundation.
8!
9!  This program is distributed in the hope that it will be useful,
10!  but WITHOUT ANY WARRANTY; without even the implied warranty of
11!  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
12!  Lesser General Public License for more details.
13!
14!  If a copy of the GNU LGPL v2.1 was not distributed with this
15!  code, you can obtain one at https://www.gnu.org/licenses/old-licenses/lgpl-2.1.en.html.
16!
17C
18C=========================================================================
19C old log:
20C010131-vebjornb: Old dalton.F split into dalton.F (program) and
21C                 dalgnr.F (subroutines).
22C951009-kenneth : Changed name from PSISTAR to DALTON
23C950621-vebjornb: First working optimization, GRAND changed to PSISTAR
24C950516-vebjornb: First effort to link all programs: GRAND
25C=========================================================================
26
27! define DALTON_VERSION
28#include "dalton_config.h"
29
30#ifdef VAR_CHEMSHELL
31      subroutine dalton
32#else
33      PROGRAM DALTON2020_DALTON
34#endif
35C
36      use parallel_models_mpi
37C
38#include "implicit.h"
39#include "dummy.h"
40#include "mtags.h"
41#include "priunit.h"
42#include "maxorb.h"
43#include "infpar.h"
44#include "infpri.h"
45#include "gnrinf.h"
46#include "inftap.h"
47#include "prpc.h"
48#include "incore.h"
49
50#if defined (VAR_MPI)
51      INCLUDE 'mpif.h'
52      LOGICAL(kind=MPI_INTEGER_KIND) :: FINISH
53      INTEGER(kind=MPI_INTEGER_KIND) :: ierr_mpi
54#endif
55#ifdef VAR_OMP
56      INTEGER omp_get_num_procs
57#endif
58      CHARACTER LABEL1*8, STHELP*10
59C defined parallel calculation types
60#include "iprtyp.h"
61
62!     this is to echo input files to the output
63      character(100)     :: line
64      logical            :: file_exists, file_open
65
66!     memory-id tags for field 1 and 2+LMWORK in work (allocated in main driver routines)
67      DATA WRKDLM/12345678D0/
68
69      ! we use this to read env variables
70      character*8 :: env_value
71
72      ! number of OMP threads
73      integer :: omp_num_threads
74C
75C
76C Start MPI, if parallel, in MPIXINIT; define MYNUM, MASTER, NODTOT
77C If not MPI, MPIXINIT defines: MYNUM  = 0; MASTER = 0; NODTOT = 0
78C
79      CALL MPIXINIT
80
81      call parallel_models_initialize_mpi()
82
83      CALL QENTER('DALTON main')
84      PARIO = .FALSE.
85C
86C     Define SLAVE variable (only false for MASTER, MYNUM.eq.0)
87C     Initialize file pointers; open error file for MASTER
88C
89      CALL GPIO_INI_DALTON
90      IF (MYNUM .EQ. MASTER) THEN
91         SLAVE = .FALSE.
92         CALL GPOPEN(LUSTAT,'DALTON.STAT','UNKNOWN',' ','FORMATTED',
93     &               IDUMMY,.FALSE.)
94      ELSE
95         SLAVE = .TRUE.
96      END IF
97
98#ifdef SYS_UNIX
99      LUERR = 0 ! stderr
100#else
101      LUERR = 6 ! stdout (which is not DALTON.OUT !)
102#endif
103
104C
105C     Read environment variables for LMWORK, LMWORK_NODES, and BASDIR
106C     If LMWORK .le. 0 use default value (LWORK_DEF)
107C
108      CALL GETMMBAS(LMWORK,BASDIR,LBASDIR)
109C
110      NBYTES  = (LMWORK+2) * 8
111      XMBYTES = (LMWORK+2) * 8.D0
112      XMBYTES = XMBYTES / (1.024D3**2)
113      GBYTES  = XMBYTES / 1.024D3
114
115      IF (MYNUM .LE. 1) THEN
116         IF (GBYTES .GE. 1.D0) THEN
117            WRITE(*,'(/A,I12,A,F7.3,A,I4)')
118     &        ' Work memory size (LMWORK+2):',
119     &        LMWORK+2,' =',GBYTES,' gigabytes; node',MYNUM
120         ELSE
121            WRITE(*,'(/A,I12,A,F8.2,A,I4)')
122     &        ' Work memory size (LMWORK+2):',
123     &        LMWORK+2,' =',XMBYTES,' megabytes; node',MYNUM
124         END IF
125         IF (MYNUM .EQ. 1 .AND. NODTOT .GT. 1) THEN
126            WRITE(*,'(A,I0,A)')
127     &        ' Same work memory size for all ',NODTOT,' slaves.'
128         END IF
129      END IF
130
131Cef begin
132C If AOSAVE .eq. .TRUE. then run the incore-version of Dalton where integrals
133C are stored in memory on each slave for reuse. Take a look at the header
134C file DALTON_SOURCE_DIR/include/incore.h  and the variables: MMWORK and IEIR.
135C MMWORK (we used MMWORK=80000000 when testing ) determines the size of
136C the chunk of array WORK which is reserved for incore-calculations.
137C The variable IEIR in incore.h must also be set to a large enough positive
138C integer (we used IEIR=3000000).
139
140      JSCORE    = 1
141      ISCORE    = JSCORE
142      LMCORE    = MMWORK
143      MMCORE    = LMCORE
144      N_SHL     = 1
145      AOSAVE    = (MMWORK.GT.1)
146      INITX     = .FALSE.
147      MSAVE     = .TRUE.
148      LINTSV    = .FALSE.
149      LINTMP    = .FALSE.
150      INDX_SHL1 = 0
151      INDX_SHL2 = 0
152      INDX_SHL3 = 0
153      INDX_SHL4 = 0
154      CALL IZERO(INDX_SHL,MXTSK)
155      IF (MMWORK .NE. MMCORE) THEN
156         CALL QUIT('MMCORE must equal MMWORK.')
157      END IF
158Cef end
159
160      LEN_BASDIR = MAX(1,LNBLNK(BASDIR))
161      IF (MYNUM .LE. 1) THEN
162         WRITE(*,'(/I4,A/5X,A)')
163     &      MYNUM,': Directories for basis set searches:',
164     &      BASDIR(1:LEN_BASDIR)
165      END IF
166C
167C     Allowing now for parallelization over different molecular geometries,
168C     we must first check if this calculation indeed is of that kind.
169C     We do this by letting the Master node read DALTON.INP, searching for
170C     the keyword .PARNMD. If found, read DALTON/MOLECULE.INP and send them
171C     to the slaves, who writes them out in their own files
172C
173c      IF (NODTOT .GT. 0) THEN
174c         IF (MYNUM .EQ. MASTER) THEN
175c            CALL PARIOT
176c         ELSE
177c            CALL PARION
178c         END IF
179c      END IF
180C
181      CALL GPOPEN(LUPRI,'DALTON.OUT','NEW',' ','FORMATTED',
182     &            IDUMMY,.FALSE.)
183      IF (SLAVE) LUSTAT = LUPRI
184C     ... GPOPEN will change name from DALTON.OUT to $WRKDIR/$OUTFIL
185C         for Master if OUTFIL is an environment variable;
186C         for slaves the name will be changed to
187C         "DALTON.OUT.n<slave no>" in GPOPEN
188
189      IF (MYNUM.EQ.MASTER) THEN
190      CALL GPOPEN(LUNDPF,'DALTON.PROP','NEW',' ','FORMATTED',
191     &            IDUMMY,.FALSE.)
192      REWIND(LUNDPF)
193      CALL GPOPEN(LUNMPF,'midasifc.prop','NEW',' ','FORMATTED',
194     &            IDUMMY,.FALSE.)
195      REWIND(LUNMPF)
196      END IF
197C
198      CALL MEMINI(LUPRI,LUERR)
199C
200      IF (MYNUM .EQ. MASTER) THEN
201C
202         CALL GPOPEN(LUCME,'DALTON.CM','NEW',' ','FORMATTED',
203     &            IDUMMY,.FALSE.)
204C
205         CALL TITLER('Dalton - An Electronic Structure Program',
206     &               '*',115)
207
208         WRITE (LUPRI,'(T5,A)') 'This is output from DALTON release '//
209     &DALTON_VERSION
210         WRITE (LUPRI,'(T10,A)')'( Web site: http://daltonprogram.org )'
211         WRITE (LUPRI,'(/3X,A/)') SEPARATOR(1:76)
212         WRITE (LUPRI,'(T5,A)')
213     &'NOTE:',
214     &' ',
215     &'Dalton is an experimental code for the evaluation of molecular',
216     &'properties using (MC)SCF, DFT, CI, and CC wave functions.',
217     &'The authors accept no responsibility for the performance of',
218     &'the code or for the correctness of the results.',
219     &' ',
220     &'The code (in whole or part) is provided under a licence and',
221     &'is not to be reproduced for further distribution without',
222     &'the written permission of the authors or their representatives.',
223     &' ',
224     &'See the home page "http://daltonprogram.org"'
225     &//' for further information.',
226     &' ',
227     &'If results obtained with this code are published,',
228     &'the appropriate citations would be both of:',
229     &' ',
230     &'   K. Aidas, C. Angeli, K. L. Bak, V. Bakken, R. Bast,',
231     &'   L. Boman, O. Christiansen, R. Cimiraglia, S. Coriani,',
232     &'   P. Dahle, E. K. Dalskov, U. Ekstroem,',
233     &'   T. Enevoldsen, J. J. Eriksen, P. Ettenhuber, B. Fernandez,',
234     &'   L. Ferrighi, H. Fliegl, L. Frediani, K. Hald, A. Halkier,',
235     &'   C. Haettig, H. Heiberg, T. Helgaker, A. C. Hennum,',
236     &'   H. Hettema, E. Hjertenaes, S. Hoest, I.-M. Hoeyvik,',
237     &'   M. F. Iozzi, B. Jansik, H. J. Aa. Jensen, D. Jonsson,',
238     &'   P. Joergensen, J. Kauczor, S. Kirpekar,',
239     &'   T. Kjaergaard, W. Klopper, S. Knecht, R. Kobayashi, H. Koch,',
240     &'   J. Kongsted, A. Krapp, K. Kristensen, A. Ligabue,',
241     &'   O. B. Lutnaes, J. I. Melo, K. V. Mikkelsen, R. H. Myhre,',
242     &'   C. Neiss, C. B. Nielsen, P. Norman, J. Olsen,',
243     &'   J. M. H. Olsen, A. Osted, M. J. Packer, F. Pawlowski,',
244     &'   T. B. Pedersen, P. F. Provasi, S. Reine, Z. Rinkevicius,',
245     &'   T. A. Ruden, K. Ruud, V. Rybkin, P. Salek, C. C. M. Samson,',
246     &'   A. Sanchez de Meras, T. Saue, S. P. A. Sauer,',
247     &'   B. Schimmelpfennig, K. Sneskov, A. H. Steindal,',
248     &'   K. O. Sylvester-Hvid, P. R. Taylor, A. M. Teale,',
249     &'   E. I. Tellgren, D. P. Tew, A. J. Thorvaldsen, L. Thoegersen,',
250     &'   O. Vahtras, M. A. Watson, D. J. D. Wilson, M. Ziolkowski',
251     &'   and H. Agren,',
252     &'   "The Dalton quantum chemistry program system",',
253     &'   WIREs Comput. Mol. Sci. 2014, 4:269–284 '//
254     &    '(doi: 10.1002/wcms.1172)',
255     &' ',
256     &'and',
257     &' ',
258     &'   Dalton, a molecular electronic structure program,',
259     &'   Release '//
260     &DALTON_VERSION
261     &//', see http://daltonprogram.org'
262C
263         WRITE (LUPRI,'(3X,A)') SEPARATOR(1:76)
264         WRITE (LUPRI,'(/T5,A//,(T3,A,T29,A,T59,A,T71,A))')
265     *'Authors in alphabetical order '//
266     &  '(Major contribution(s) in parenthesis):',
267     *'Kestutis Aidas,','Vilnius University,','Lithuania',
268     &  '(QM/MM)',
269     *'Celestino Angeli,','University of Ferrara,','Italy',
270     &  '(NEVPT2)',
271     *'Keld L. Bak,','UNI-C,','Denmark',
272     &  '(AOSOPPA, non-adiabatic coupling, magnetic properties)',
273     *'Vebjoern Bakken,','University of Oslo,','Norway',
274     &  '(DALTON; geometry optimizer, symmetry detection)',
275     *'Radovan Bast,','UiT The Arctic U. of Norway,','Norway',
276     &  '(DALTON installation and execution frameworks)',
277     *'Pablo Baudin,','University of Valencia,','Spain',
278     &  '(Cholesky excitation energies)',
279     *'Linus Boman,','NTNU,', 'Norway',
280     &  '(Cholesky decomposition and subsystems)',
281     *'Bruno N. Cabral Tenorio,','Technical Univ. of Denmark,',
282     &'Denmark', '(Lanczos-based CC methods)',
283     *'Ove Christiansen,','Aarhus University,','Denmark',
284     &  '(CC module)',
285     *'Renzo Cimiraglia,','University of Ferrara,','Italy',
286     &  '(NEVPT2)',
287     *'Sonia Coriani,','Technical Univ. of Denmark,','Denmark',
288     &  '(CC and RESPONS modules)',
289     *'Janusz Cukras,','University of Warsaw,','Poland',
290     &  '(MChD in RESPONS)',
291     *'Paal Dahle,','University of Oslo,','Norway',
292     &  '(Parallelization)',
293     *'Erik K. Dalskov,','UNI-C,','Denmark',
294     &  '(SOPPA)',
295     *'Thomas Enevoldsen,','Univ. of Southern Denmark,','Denmark',
296     &  '(SOPPA)',
297     *'Janus J. Eriksen,','Aarhus University,','Denmark',
298     &  '(Polarizable embedding, TDA)',
299     *'Rasmus Faber,','Technical Univ. of Denmark,','Denmark',
300     &  '(Vib.avg. NMR w. SOPPA, parallel AO-SOPPA, CVS-CC and EOMCC)',
301     *'Tobias Fahleson,','KTH Stockholm,','Sweden',
302     &  '(Damped cubic response)',
303     *'Berta Fernandez,','U. of Santiago de Compostela,','Spain',
304     &  '(Doublet spin, ESR in RESPONS)',
305     *'Lara Ferrighi,','Aarhus University,','Denmark',
306     &  '(PCM Cubic response)',
307     *'Heike Fliegl,','University of Oslo,','Norway',
308     &  '(CCSD(R12))',
309     *'Luca Frediani,','UiT The Arctic U. of Norway,','Norway',
310     &  '(PCM)',
311     *'Emmanuel Fromager,','Univ. Strasbourg,','France',
312     &  '(MC-srDFT)',
313     *'Bin Gao,','UiT The Arctic U. of Norway,','Norway',
314     &  '(Gen1Int library)',
315     *'Andre S. P. Gomes,','CNRS/Universite de Lille,','France',
316     &  '(Frozen density embedding)',
317     *'Christof Haettig,','Ruhr-University Bochum,','Germany',
318     &  '(CC module)',
319     *'Kasper Hald,', 'Aarhus University,','Denmark',
320     &  '(CC module)',
321     *'Asger Halkier,','Aarhus University,','Denmark',
322     &  '(CC module)',
323     *'Frederik Beyer Hansen,','University of Copenhagen,','Denmark',
324     &  '(Parallel AO-SOPPA)',
325     *'Erik D. Hedegaard,','Univ. of Southern Denmark,','Denmark',
326     & '(Polarizable embedding, QM/MM, and MC-srDFT)',
327     *'Hanne Heiberg,','University of Oslo,','Norway',
328     &  '(Geometry analysis, selected one-electron integrals)',
329     *'Trygve Helgaker,','University of Oslo,','Norway',
330     &  '(DALTON; ABACUS, ERI, DFT modules, London, and much more)',
331     *'Alf Christian Hennum,','University of Oslo,','Norway',
332     &  '(Parity violation)',
333     *'Hinne Hettema,','University of Auckland,','New Zealand',
334     &  '(Quadratic response in RESPONS; SIRIUS supersymmetry)',
335     *'Eirik Hjertenaes,', 'NTNU,','Norway',
336     &  '(Cholesky decomposition)',
337     *'Pi A. B. Haase,','University of Copenhagen,','Denmark',
338     &  '(Triplet AO-SOPPA)',
339     *'Maria Francesca Iozzi,','University of Oslo,','Norway',
340     &  '(RPA)',
341     *'Christoph Jacob','TU Braunschweig','Germany',
342     &  '(Frozen density embedding)',
343     *'Brano Jansik','Technical Univ. of Ostrava','Czech Rep.',
344     &  '(DFT cubic response)',
345     *'Hans Joergen Aa. Jensen,','Univ. of Southern Denmark,','Denmark',
346     &  '(DALTON; SIRIUS, RESPONS, ABACUS, MC-srDFT, and much more)'
347     *,'Dan Jonsson,','UiT The Arctic U. of Norway,','Norway',
348     &  '(Cubic response in RESPONS module)',
349     *'Poul Joergensen,','Aarhus University,','Denmark',
350     &  '(RESPONS, ABACUS, and CC modules)',
351     *'Maciej Kaminski,','University of Warsaw,','Poland',
352     &  '(CPPh in RESPONS)',
353     *'Joanna Kauczor,','Linkoeping University,','Sweden',
354     &  '(Complex polarization propagator (CPP) module)',
355     *'Sheela Kirpekar,','Univ. of Southern Denmark,','Denmark',
356     &  '(Mass-velocity & Darwin integrals)',
357     *'Erik Kjellgren,','Univ. of Southern Denmark,','Denmark',
358     &  '(MC-srDFT)',
359     *'Wim Klopper,','KIT Karlsruhe,','Germany',
360     &  '(R12 code in CC, SIRIUS, and ABACUS modules)',
361     *'Stefan Knecht,','ETH Zurich,','Switzerland',
362     &  '(Parallel CI and MCSCF, MC-srDFT)',
363     *'Rika Kobayashi,','Australian National Univ.,','Australia',
364     &  '(DIIS in CC, London in MCSCF)',
365     *'Henrik Koch,','NTNU,', 'Norway',
366     &  '(CC module, Cholesky decomposition)',
367     *'Jacob Kongsted,','Univ. of Southern Denmark,','Denmark',
368     &  '(Polarizable embedding, polarizable density embedding, QM/MM)',
369     *'Andrea Ligabue,','University of Modena,','Italy',
370     &  '(CTOCD, AOSOPPA)',
371     *'Nanna H. List','Univ. of Southern Denmark,','Denmark',
372     &  '(Polarizable embedding)',
373     *'Ola B. Lutnaes,','University of Oslo,','Norway',
374     &  '(DFT Hessian)',
375     *'Juan I. Melo,','University of Buenos Aires,','Argentina',
376     &  '(LRESC, Relativistic Effects on NMR Shieldings)',
377     *'Kurt V. Mikkelsen,','University of Copenhagen,','Denmark',
378     &  '(MC-SCRF and QM/MM)',
379     *'Rolf H. Myhre,', 'NTNU,','Norway',
380     &  '(Subsystems and CC3)',
381     *'Christian Neiss,','Univ. Erlangen-Nuernberg,','Germany',
382     &  '(CCSD(R12))',
383     *'Christian B. Nielsen,','University of Copenhagen,','Denmark',
384     &  '(QM/MM)',
385     *'Patrick Norman,','KTH Stockholm,','Sweden',
386     &  '(Cubic response and complex frequency response in RESPONS)',
387     *'Jeppe Olsen,','Aarhus University,','Denmark',
388     &  '(SIRIUS CI/density modules)',
389     *'Jogvan Magnus H. Olsen,','Aarhus University,','Denmark',
390     &  '(Polarizable embedding and polarizable density embedding)',
391     *'Anders Osted,','Copenhagen University,','Denmark',
392     &  '(QM/MM)',
393     *'Martin J. Packer,','University of Sheffield,','UK',
394     &  '(SOPPA)',
395     *'Filip Pawlowski,','Kazimierz Wielki University,','Poland',
396     &  '(CC3)',
397     *'Morten N. Pedersen,','Univ. of Southern Denmark,','Denmark',
398     &  '(Polarizable embedding)',
399     *'Thomas B. Pedersen,','University of Oslo,','Norway',
400     &  '(Cholesky decomposition)',
401     *'Patricio F. Provasi,','University of Northeastern,','Argentina',
402     &  '(Analysis of coupling constants in localized orbitals)',
403     *'Peter Reinholdt,','Univ. of Southern Denmark,','Denmark',
404     &  '(Polarizable embedding and polarizable density embedding)',
405     *'Zilvinas Rinkevicius,','KTH Stockholm,','Sweden',
406     &  '(open-shell DFT, ESR)',
407     *'Elias Rudberg,','KTH Stockholm,','Sweden',
408     &  '(DFT grid and basis info)',
409     *'Torgeir A. Ruden,','University of Oslo,','Norway',
410     &  '(Numerical derivatives in ABACUS)',
411     *'Kenneth Ruud,','UiT The Arctic U. of Norway,','Norway',
412     &  '(DALTON; ABACUS magnetic properties and much more)',
413     *'Pawel Salek,','KTH Stockholm,','Sweden',
414     &  '(DALTON; DFT code)',
415     *'Claire C. M. Samson','University of Karlsruhe','Germany',
416     &  '(Boys localization, r12 integrals in ERI)',
417     *'Alfredo Sanchez de Meras,','University of Valencia,','Spain',
418     &  '(CC module, Cholesky decomposition)',
419     *'Trond Saue,','Paul Sabatier University,','France',
420     &  '(direct Fock matrix construction)',
421     *'Stephan P. A. Sauer,','University of Copenhagen,','Denmark',
422     &  '(SOPPA(CCSD), SOPPA prop., AOSOPPA, vibrational g-factors)',
423     *'Bernd Schimmelpfennig,','Forschungszentrum Karlsruhe,',
424     &  'Germany','(AMFI module)',
425     *'Anna K. Schnack-Petersen,','University of Copenhagen,','Denmark',
426     &  '(RPA(D) and HRPA(D) properties in AOSOPPA)',
427     *'Kristian Sneskov,','Aarhus University,','Denmark',
428     &  '(Polarizable embedding, QM/MM)',
429     *'Arnfinn H. Steindal,','UiT The Arctic U. of Norway,',
430     &  'Norway','(parallel QM/MM, polarizable embedding)',
431     *'Casper Steinmann,','Aalborg University,','Denmark',
432     &  '(QFIT, polarizable embedding)',
433     *'K. O. Sylvester-Hvid,','University of Copenhagen,','Denmark',
434     &  '(MC-SCRF)',
435     *'Peter R. Taylor,','VLSCI/Univ. of Melbourne,','Australia',
436     &  '(Symmetry handling ABACUS, integral transformation)',
437     *'Andrew M. Teale,','University of Nottingham,','England',
438     &  '(DFT-AC, DFT-D)',
439     *'David P. Tew,','University of Bristol,','England',
440     &  '(CCSD(R12))',
441     *'Julien Toulouse,','Sorbonne Univ.,','France',
442     &  '(MC-srDFT)',
443     *'Olav Vahtras,','KTH Stockholm,','Sweden',
444     &  '(Triplet response, spin-orbit, ESR, TDDFT, open-shell DFT)',
445     *'Lucas Visscher,','Vrije Universiteit Amsterdam,','Netherlands',
446     &  '(Frozen density embedding)',
447     *'David J. Wilson,','La Trobe University,','Australia',
448     &  '(DFT Hessian and DFT magnetizabilities)',
449     *'Hans Agren,','KTH Stockholm,','Sweden',
450     &  '(SIRIUS module, RESPONS, MC-SCRF solvation)'
451C
452         WRITE (LUPRI,'(1X,A)') SEPARATOR
453C
454C     Stamp date and time and hostname to output
455C
456         CALL TSTAMP('INIT',LUPRI)
457         IF (GBYTES .GE. 1.D0) THEN
458           WRITE(LUPRI,'(/A,I12,A,F7.3,A)')
459     &        ' * Work memory size             :',
460     &        LMWORK,' =',GBYTES,' gigabytes.'
461         ELSE
462           WRITE(LUPRI,'(/A,I12,A,F8.2,A)')
463     &        ' * Work memory size             :',
464     &        LMWORK,' =',XMBYTES,' megabytes.'
465         END IF
466         IF (AOSAVE) WRITE(LUPRI,'(A,I12)')
467     &        ' + memory for in-core integrals :',MMCORE
468C
469         WRITE(LUPRI,'(/A)')
470     &   ' * Directories for basis set searches:'
471         J_END = LNBLNK(BASDIR)
472         IF (J_END .LE. 0) THEN
473            WRITE(LUPRI,'(5X,A)') 'NONE SPECIFIED!'
474         ELSE
475C           structure of BASDIR: 'dir1:dir2:dir3'
476            I = 0
477            I_END = -1
478   20       CONTINUE
479               I_ST = I_END + 2
480               I_END = INDEX(BASDIR(I_ST:J_END),':') - 1
481               IF (I_END .LE. 0) THEN
482                  I_END = J_END
483               ELSE
484                  I_END = I_ST - 1 + I_END
485               END IF
486               IF (I_END .GT. I_ST) THEN
487                  I = I + 1
488                  WRITE(LUPRI,'(I4,2A)') I,') ',BASDIR(I_ST:I_END)
489                  GO TO 20
490               END IF
491         END IF
492
493         write(lupri, '(//a/a/)')
494     &     'Compilation information',
495     &     '-----------------------'
496         call print_binary_info(lupri)
497
498#if defined (VAR_MPI)
499         IF (NODTOT .GT. 0) THEN
500            WRITE(LUPRI,'(/A,I0,A)')
501     &      " * MPI parallel run using ", NODTOT+1, " processes."
502         ELSE
503            WRITE(LUPRI,'(/A,I5,A)')
504     &      " * Sequential calculation."
505         END IF
506#endif
507
508#ifdef VAR_OMP
509         call getenv('OMP_NUM_THREADS', env_value)
510         read(env_value, '(i6)') omp_num_threads
511         if (omp_num_threads > 1) then
512            WRITE(LUPRI,'(/A,I0,A)')
513     &      " * OpenMP run using ", omp_num_threads, " threads."
514         end if
515#endif
516
517!        print content of input file to output
518         inquire(file = 'DALTON.INP', exist = file_exists,
519     &           opened = file_open)
520         if (file_exists) then
521            if (.not. file_open) then
522               LUCMD = -1
523               CALL GPOPEN(LUCMD,'DALTON.INP','OLD',' ','FORMATTED',
524     &                     IDUMMY,.FALSE.)
525            end if
526            rewind(LUCMD)
527         else
528            call quit('DALTON.INP file not found')
529         end if
530         write(lupri, '(//A/A/)')
531     &      '   Content of the .dal input file',
532     &      ' ----------------------------------'
533         do while (.true.)
534            read(LUCMD, '(a)', end=31) line
535            write(lupri, '(a)') trim(line)
536         end do
537 31      continue
538         if (.not. file_open) CALL GPCLOSE(LUCMD,'KEEP')
539
540!        print content of molecule file to output
541         inquire(file = 'MOLECULE.INP', exist = file_exists)
542         if (file_exists) then
543            LUMOL = -1
544            CALL GPOPEN(LUMOL,'MOLECULE.INP','OLD',' ','FORMATTED',
545     &                  IDUMMY,.FALSE.)
546            write(lupri, '(//a/a/)')
547     &         '   Content of the .mol file',
548     &         ' ----------------------------'
549            rewind(LUMOL)
550            do while (.true.)
551               read(LUMOL, '(a)', end=32) line
552               write(lupri, '(a)') trim(line)
553            end do
554 32         continue
555            CALL GPCLOSE(LUMOL,'KEEP')
556         else
557!           call quit('MOLECUL.INP file not found')
558! somebody forgot with this 'call quit' that it is OK
559! to have the "MOLECULE.INP" input inside the "DALTON.INP" file,
560! thus the *nevpt2* tests fail with this 'call quit' uncommented.
561         end if
562
563C        Read general input for Dalton
564C
565         CALL DALTON_GNRLINP
566
567#if defined (VAR_MPI)
568         CALL RECVNAMES
569
570         CALL MOLDEN_HEAD
571!
572!        call main driver routine
573         call dalton_exedrv(lmwork,wrkdlm)
574
575         CALL MOLDEN_TAIL
576
577         CALL FLSHFO(LUPRI)
578
579         LABEL1 = 'ALL_DONE'
580         STHELP = 'THE_END   '
581         ZERO   = 0.0D0
582         CALL WRIPRO(0.0D0,STHELP,666,
583     *               LABEL1,LABEL1,LABEL1,LABEL1,
584     *               ZERO,ZERO,ZERO,1,0,0,0)
585         IPRLOC = MAX(IPRUSR,IPRRHF)
586C        IF (IPRLOC .GT. 2) CALL PRPRPC(LUNDPF,1,DUMMY,NPRMI)
587         IF (IPRLOC .GT. 2) CALL PRPRPC(LUNMPF,1,DUMMY,NPRMI)
588C
589C
590C     Send ending signal to slaves
591C
592         IPRTYP = ENDING_SIGNAL
593         CALL MPIXBCAST(IPRTYP,1,'INTEGER',MYNUM)
594
595         CALL GPCLOSE(LUCME,'KEEP')
596         CALL GPCLOSE(LUNDPF,'KEEP')
597         CALL GPCLOSE(LUNMPF,'KEEP')
598         CALL GPCLOSE(LUPRI,'KEEP')
599         CALL GPCLOSE(LUSTAT,'KEEP')
600      ELSE
601C     ... I am a "slave" ... (MYNUM .ne. MASTER)
602C ahs only print if number of cores are less than 16
603!        IF (NODTOT .LE. 15) THEN
604!        hjaaj May 2017: instead print from max 15 slaves
605!           spaced evenly out among all slaves 1..NODTOT
606!        ITEST = (NODTOT-1)/15 + 1
607!        IF (MYNUM .EQ. 1 .OR. MOD(MYNUM,ITEST) .EQ. 0) THEN
608           CALL TITLER('DALTON - An electronic structure program',
609     &               '*',111)
610           WRITE(LUPRI,'(/A)') ' * Release: '//
611     &     DALTON_VERSION
612           WRITE(LUPRI,'(/A,I0,A,I0)')
613     &     ' * Output from slave no. ',MYNUM,' out of ',NODTOT
614           CALL TSTAMP('INIT',LUPRI)
615           IF (GBYTES .GE. 1.D0) THEN
616             WRITE(LUPRI,'(/A,I12,A,F7.3,A)')' * Work memory size :',
617     &        LMWORK,' =',GBYTES,' gigabytes.'
618           ELSE
619             WRITE(LUPRI,'(/A,I12,A,F8.2,A)')' * Work memory size :',
620     &        LMWORK,' =',XMBYTES,' megabytes.'
621           END IF
622           IF (AOSAVE) THEN
623             WRITE(LUPRI,'(A,I12)')' + memory for in-core integrals :',
624     &        MMCORE
625           END IF
626           CALL FLSHFO(LUPRI)
627!        END IF
628         IF (PARIO) THEN
629            CALL DALTON_GNRLINP
630
631            CALL MOLDEN_HEAD
632
633!           call main driver routine (co-workers do it as well if PARIO == .true.)
634            call dalton_exedrv(lmwork,wrkdlm)
635
636            CALL MOLDEN_TAIL
637
638            LABEL1 = 'ALL_DONE'
639            STHELP = 'THE_END   '
640            ZERO   = 0.0D0
641            CALL WRIPRO(0.0D0,STHELP,666,
642     *                  LABEL1,LABEL1,LABEL1,LABEL1,
643     *                  ZERO,ZERO,ZERO,1,0,0,0)
644            IPRLOC = MAX(IPRUSR,IPRRHF)
645C           IF (IPRLOC .GT. 2) CALL PRPRPC(LUNDPF,1,DUMMY,NPRMI)
646            IF (IPRLOC .GT. 2) CALL PRPRPC(LUNMPF,1,DUMMY,NPRMI)
647C
648C     Send ending signal to slaves
649C
650            FINISH = .TRUE.
651            count_mpi = 1
652            master_mpi = MASTER
653            CALL MPI_BCAST(FINISH,count_mpi,my_MPI_LOGICAL,master_mpi,
654     &           MPI_COMM_WORLD,err_mpi)
655            CALL GPCLOSE(LUCME,'KEEP')
656            CALL GPCLOSE(LUNDPF,'KEEP')
657            CALL GPCLOSE(LUNMPF,'KEEP')
658            CALL GPCLOSE(LUSTAT,'KEEP')
659            LUSTAT = LUPRI ! LUSTAT may be used in GPCLOSE
660            CALL GPCLOSE(LUPRI,'KEEP')
661         ELSE
662            CALL SENDNAME(MASTER,MYNUM)
663
664!           call main co-workers driver routine
665            call dalton_nodedriver(lmwork,wrkdlm)
666         END IF
667      END IF
668
669#if defined(BUILD_GEN1INT)
670C...  added by Bin Gao, May 13, 2012
671C...  terminates Gen1Int interface
672      call gen1int_host_finalize()
673#endif
674
675      call parallel_models_finalize_mpi(nodtot+1)
676
677      CALL MPIXFINALIZE
678C
679C     All below is not relevant for MPI as all proc. has stopped
680C     with CALL MPI_FINALIZE(err_mpi)
681C
682#else /* VAR_MPI */
683!     sequential build
684      END IF ! mynum .eq. master
685
686      CALL MOLDEN_HEAD
687
688!     call main co-workers driver routine
689      call dalton_exedrv(lmwork,wrkdlm)
690
691      CALL MOLDEN_TAIL
692
693      LABEL1 = 'ALL_DONE'
694      STHELP = 'THE_END   '
695      ZERO   = 0.0D0
696      CALL WRIPRO(0.0D0,STHELP,666,
697     *            LABEL1,LABEL1,LABEL1,LABEL1,
698     *            ZERO,ZERO,ZERO,1,0,0,0)
699      IPRLOC = MAX(IPRUSR,IPRRHF)
700C     IF (IPRLOC .GT. 2) CALL PRPRPC(LUNDPF,1,DUMMY,NPRMI)
701      IF (IPRLOC .GT. 2) CALL PRPRPC(LUNMPF,1,DUMMY,NPRMI)
702      CALL GPCLOSE(LUCME,'KEEP')
703      CALL GPCLOSE(LUNDPF,'KEEP')
704      CALL GPCLOSE(LUNMPF,'KEEP')
705      CALL GPCLOSE(LUSTAT,'KEEP')
706      LUSTAT = LUPRI ! LUSTAT may be used in GPCLOSE
707      CALL GPCLOSE(LUPRI,'KEEP')
708#if defined(BUILD_GEN1INT)
709C...  added by Bin Gao, Oct. 2, 2011
710C...  terminates Gen1Int interface
711      call gen1int_host_finalize()
712#endif
713#endif
714
715
716      ! radovan: close (hopefully) all open files
717      !          otherwise we run out of units
718      !          if we execute Dalton many times in sequence
719      !          (ChemShell for example)
720      do iunit = 1, 99
721         inquire(iunit, opened=file_open)
722         if (file_open) close(iunit, status='keep')
723      end do
724      ! clear the corresponding table and hope for the best
725      call clear_iuntab()
726
727      CALL QEXIT('DALTON main')
728      END
729
730      SUBROUTINE GETMMBAS(LMWORK,BASDIR,LBASDIR)
731
732!     module dependencies
733      use parallel_communication_models_mpi
734
735#include "implicit.h"
736#include "priunit.h"
737      CHARACTER*(*) BASDIR
738      CHARACTER*20  WRKMEM, NODE_WRKMEM
739#include "maxorb.h"
740#include "infpar.h"
741!     this is to scan the input file for asynchronous memory allocation (RMA model parallelization)
742      character(100)     :: line
743      logical            :: file_exists, file_open, rma_allocation
744C     If LMWORK .le. 0 use default value (LWORK_DEF)
745      PARAMETER (LWORK_DEF = INSTALL_WRKMEM)
746
747      rma_allocation = .false.
748C
749C     Dynamic memory allocation; check for WRKMEM and NODE_WRKMEM in environment
750C
751      IF (.NOT. SLAVE) THEN
752C     ... only master reads environment, transferred
753C         to slave with MPI_BCAST.
754C         On some MPI implementations it requires extra
755C         options to transfer UNIX environment variables to
756C         slave nodes. /hjaaj June 2005
757        WRKMEM = ' '
758        NODE_WRKMEM = ' '
759        CALL GETENV('WRKMEM',WRKMEM)
760        CALL GETENV('NODE_WRKMEM',NODE_WRKMEM)
761        CALL GETENV('BASDIR',BASDIR)
762        !debug! write (0,*) 'GETENV WRKMEM ',WRKMEM
763        !debug! write (0,*) 'GETENV BASDIR ',BASDIR
764
765!       WRKMEM -> LMWORK
766
767        READ(WRKMEM, '(I20)', ERR=10) LMWORK
768        IF (LMWORK .GE. 0) GO TO 20
769   10     WRITE(LUERR,'(/3A/A,I20)')
770     &    ' DALTON: WRKMEM conversion error; WRKMEM = "',WRKMEM,'"',
771     &    ' DALTON: read as LMWORK = ',LMWORK
772          LMWORK = 0
773   20   CONTINUE
774C
775C       Either default memory size,
776C
777        IF (LMWORK .LE. 0) THEN
778          LMWORK = LWORK_DEF
779          WRITE (*,'(/A,I20)')
780     &      ' DALTON: default work memory size used.',LMWORK
781        ELSE
782C
783C       ... or user specified memory size.
784C
785          LWMEM = MAX(1,LNBLNK(WRKMEM))
786          WRITE(*,'(/A/3A)')
787     &    ' DALTON: user specified work memory size used,',
788     &    '         environment variable WRKMEM = "',WRKMEM(1:LWMEM),'"'
789        END IF
790
791!       NODE_WRKMEM -> LMWORK_NODES
792
793      IF (NODTOT .GT. 0) THEN ! not relevant if no slaves
794        READ(NODE_WRKMEM, '(I20)', ERR=60) LMWORK_NODES
795        IF (LMWORK_NODES .GE. 0) GO TO 70
796   60     WRITE(LUERR,'(/3A/A,I20)')
797     &    ' DALTON: NODE_WRKMEM conversion error; NODE_WRKMEM = "',
798     &      NODE_WRKMEM,'"',
799     &    ' DALTON: read as LMWORK_NODES = ',LMWORK_NODES
800          LMWORK_NODES = 0
801   70   CONTINUE
802C
803C       Either default memory size,
804C
805        IF (LMWORK_NODES .LE. 0) THEN
806          LMWORK_NODES = LMWORK
807          WRITE (*,'(/A,I20)')
808     &    ' DALTON: master work memory size also used for slaves.',
809     &    LMWORK_NODES
810        ELSE
811C
812C       ... or user specified memory size.
813C
814          LWMEM = MAX(1,LNBLNK(NODE_WRKMEM))
815          WRITE(*,'(/A/3A)')
816     &    ' DALTON: user specified work memory size used,',
817     &    '         environment variable NODE_WRKMEM = "',
818     &    NODE_WRKMEM(1:LWMEM),'"'
819        END IF
820      END IF ! NODTOT .gt. 0
821
822      END IF ! .not. slave
823
824      IF (LNBLNK(BASDIR) .LT. 1) BASDIR = './:'
825#ifdef VAR_MPI
826
827!     scan input for RMA (remote memory access model parallelization) keyword
828!     reason: we need asynchronous memory allocation as this is our default RMA scheme in:
829!             --> hermit
830!             --> mcscf (soon)
831!             --> lucita/gasci (soon)
832      if(mynum == master)then
833!       print content of input file to output
834        inquire(file='DALTON.INP',exist=file_exists,opened=file_open)
835        if(file_exists) then
836          if(.not. file_open) then
837            LUCMD = -1
838            CALL GPOPEN(LUCMD,'DALTON.INP','OLD',' ','FORMATTED',
839     &                  IDUMMY,.FALSE.)
840          end if
841            rewind(LUCMD)
842         else
843            call quit('GETMMBAS: DALTON.INP file not found')
844        end if
845!       scan for .RMA-MD
846        do while (.true.)
847           read(LUCMD, '(a50)', end=81) line
848           if(line(1:7) == '.RMA-MD')then
849             rma_allocation = .true.
850             goto 81
851           end if
852        end do
853 81     continue
854        if (.not. file_open) CALL GPCLOSE(LUCMD,'KEEP')
855      end if
856!     end of keyword scan
857
858!     Transfer LMWORK_NODES and BASDIR to slaves
859!     ------------------------------------------
860
861      CALL MPIXBCAST(LMWORK_NODES,1,'INTEGE',MASTER)
862      CALL MPIXBCAST(LMWORK      ,1,'INTEGE',MASTER)
863      CALL MPIXBCAST(rma_allocation,1,'LOGICA',MASTER)
864      CALL MPIXBCAST(BASDIR,LBASDIR,'STRING',MASTER)
865
866      IF(SLAVE)THEN
867        if(rma_allocation)then
868!         node/NUMA master get the same memory as the master only real slaves get less
869          if(communication_info_mpi%my_shmem_node_id /= 0)then
870            LMWORK = LMWORK_NODES
871          end if
872        else
873          LMWORK = LMWORK_NODES
874        end if
875      END IF
876#endif /* defined(VAR_MPI) */
877      END
878C  /* Deck GPIO_INI_DALTON */
879      SUBROUTINE GPIO_INI_DALTON
880C
881C 1-Mar-2000 K.Ruud
882C 10-Apr-2000 HJAaJ: changed initialization to -9xxx from 0,
883C                    because unit 0 is stderr on many systems.
884C                    The number -9xxx is used for easier debugging.
885C
886#include "mxcent.h"
887#include "eritap.h"
888#include "priunit.h"
889#include "inftap.h"
890#include "ccinftap.h"
891#include "r12int.h"
892C
893C     priunit.h:
894C
895      LUERR  = -8000
896      LUSTAT = -8001
897      LUW4   = -8004
898      LUCMD  = -8005
899      LUPRI  = -8006
900      LUPOT  = -8009
901      NINFO  = 0
902      NWARN  = 0
903      IPRSTAT= 0
904C
905C     inftap.h:
906C
907      LUCME  = -9003
908      LUMOL  = -9004
909      LUPROP = -9005
910      LUSOL  = -9006
911      LUINTA = -9007
912      LUONEL = -9008
913      LUSUPM = -9009
914      LUTLM  = -9010
915      LUDA1  = -9011
916      LUITMP = -9012
917      LU2DER = -9013
918      LUDASP = -9014
919      LURDR  = -9016
920      LURDI  = -9017
921      LUGDR  = -9018
922      LUGDI  = -9019
923      LUGDT  = -9020
924      LURDT  = -9021
925      LUDFCK = -9022
926      LUSFDA = -9023
927      LUFDC  = -9024
928      LUWLK  = -9025
929      LUPAO  = -9026
930      LUPAS  = -9027
931      LUNR1  = -9028
932      LUNR3  = -9029
933      LUNR5  = -9030
934      LUINTR = -9032
935      LUPMOM = -9034
936      LUMOM  = -9035
937      LUEIND = -9036
938      LUENUC = -9037
939      LUESITE= -9038
940      LUEOBAR= -9039
941      LUVDWSE= -9040
942      LUENSA = -9041
943      LUQM3E = -9042
944      LUQM3P = -9043
945      LUOSCR = -9044
946      LUMMOV = -9045
947      LUOVER = -9046
948      LUNDPF = -9047
949      LUSRINT =-9048
950      LUNMPF = -9049
951      LUMOLDEN    =-9050
952      LUMOLDEN_MOS=-9051
953C
954C     /R12INT/ (WK/UniKA/26-11-2002).
955C
956      DO I = 1, 5
957         LUR12(I) = -9100 - I
958      ENDDO
959C
960C     /ERITAP/
961C
962      DO I = 0, MXCOOR
963         LUAORC(I) = -10000 - I
964      ENDDO
965C
966C     /RSPTAP/
967C
968      LUAHSO = -9201
969      LUCRV1 = -9202
970      LUCRV2 = -9203
971      LUXYVE = -9204
972      LUCRVE = -9205
973      LURSP3 = -9206
974      LURSP4 = -9207
975      LURSP5 = -9208
976      LUMHSO = -9209
977      LURSP  = -9210
978C
979C     /SIRTAP/
980C
981      LUINTM = -9301
982      LUIT1  = -9302
983      LUIT2  = -9303
984      LUIT3  = -9304
985      LUIT5  = -9305
986      LUINF  = -9306
987      LUH2AC = -9307
988      LUSIFC = -9308
989C     next two for NEWTRA
990      LUORDA = -9309
991      LUMINT = -9310
992C
993      FNSOL  = 'AOSOLINT'
994      ABARDR = 'ABACUS.RD'
995      ABARDI = 'ABACUS.RDI'
996      ABAGDR = 'ABACUS.GD'
997      ABAGDI = 'ABACUS.GDI'
998      ABAGDT = 'ABACUS.GDT'
999      ABARDT = 'ABACUS.RDT'
1000      ABADFK = 'ABACUS.DFK'
1001      ABASF  = 'ABACUS.SF'
1002      ABATLM = 'ABACUS.TLM'
1003      ABAWLK = 'DALTON.WLK'
1004      ABAIRC = 'DALTON.IRC'
1005      ABATRJ = 'DALTON.TRJ'
1006      ABANR1 = 'ABAENR.RST'
1007      ABANR3 = 'ABAENR.BVC'
1008      ABANR5 = 'ABAENR.SVC'
1009      FNINTM = 'MOTWOINT'
1010      FNSUPM = 'AOSUPINT'
1011      FNONEL = 'AOONEINT'
1012      FNSIFC = 'SIRIFC'
1013      LBSIFC = 'SIR IPH '
1014C
1015      RETURN
1016      END
1017