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!
17!
18C
19C  /* Deck hypinp */
20      SUBROUTINE HYPINP(WORD)
21C
22      use pelib_interface, only: use_pelib
23#include "implicit.h"
24C
25#include "priunit.h"
26#include "infrsp.h"
27#include "rspprp.h"
28#include "inflr.h"
29#include "infhyp.h"
30#include "infpri.h"
31#include "infspi.h"
32#include "maxorb.h"
33#include "infinp.h"
34#include "inftap.h"
35#include "infrank.h"
36#include "inforb.h"
37! gnrinf.h: QM3, THR_REDFAC
38#include "gnrinf.h"
39C
40      LOGICAL NEWDEF, BCFREQ, COMFREQ, NXTLAB
41      PARAMETER ( NTABLE = 30, ZERO = 0.0D0 , THRFRQ = 1.0D-14 )
42      CHARACTER PROMPT*1, WORD*7, TABLE(NTABLE)*7, WORD1*7
43      CHARACTER*8 LABEL, LABFND, RTNLBL(2)
44C
45      DATA TABLE /'.BPROP ', '.BFREQ ', '.CPROP ', '.CFREQ ',
46     *            '.APROP ', '.MAX IT', '.THCLR ', '.MAXITO',
47     *            '.PRINT ', '.DIPLEN', '.E3TEST', '.XXXXXX',
48     *            '.ISPABC', '.A2TEST', '.X2TEST', '.REFCHK',
49     *            '.TSTJEP', '.SHG   ', '.POCKEL', '.OPTREF',
50     *            '.FREQUE', '.DIPLNX', '.DIPLNY', '.DIPLNZ',
51     *            '.ASPIN ', '.BSPIN ', '.CSPIN ', '.SOSHIE',
52     *            '.SOSPIN', '.QLOP '/
53C
54C READ quadratic response INPUT
55C
56      NEWDEF = (WORD .EQ. '*QUADRA')
57      ICHANG = 0
58      IF (NEWDEF) THEN
59      HYPCAL = .TRUE.
60      BCFREQ = .FALSE.
61      COMFREQ= .FALSE.
62      QLOP = .FALSE.
63      WORD1 = WORD
64  100    CONTINUE
65            READ (LUCMD, '(A7)') WORD
66            CALL UPCASE(WORD)
67            PROMPT = WORD(1:1)
68            IF (PROMPT .EQ. '!' .OR. PROMPT .EQ. '#') GO TO 100
69            IF (PROMPT .EQ. '.') THEN
70               ICHANG = ICHANG + 1
71               DO 200 I = 1, NTABLE
72                  IF (TABLE(I) .EQ. WORD) THEN
73                     GO TO (1,2,3,4,5,6,7,8,9,10,11,12,13,
74     *                    14,15,16,17,18,19,20,21,22,23,24,
75     *                    25,26,27,28,29,30),I
76                  END IF
77  200          CONTINUE
78               IF (WORD .EQ. '.OPTION') THEN
79                 CALL PRTAB(NTABLE,TABLE,WORD1//' input keywords',LUPRI)
80                 GO TO 100
81               END IF
82               WRITE (LUPRI,'(/,3A,/)') ' KEYWORD "',WORD,
83     *            '" NOT RECOGNIZED IN RSPQR.'
84               CALL PRTAB(NTABLE,TABLE,WORD1//' input keywords',LUPRI)
85               CALL QUIT(' ILLEGAL KEYWORD IN *QUADRA')
86    1          CONTINUE
87                  READ(LUCMD,'(BN,A,I8)')LABEL, IRANKB
88                  IF (LABEL .EQ. 'ANGMOM  ') THEN
89                     BQROP(INDPRP('XANGMOM ')) = .TRUE.
90                     BQROP(INDPRP('YANGMOM ')) = .TRUE.
91                     BQROP(INDPRP('ZANGMOM ')) = .TRUE.
92                  ELSE IF (LABEL .EQ. '1SPNORB ') THEN
93                     BQROP(INDPRP('X1SPNORB')) = .TRUE.
94                     BQROP(INDPRP('Y1SPNORB')) = .TRUE.
95                     BQROP(INDPRP('Z1SPNORB')) = .TRUE.
96                  ELSE IF (LABEL .EQ. '2SPNORB ') THEN
97                     BQROP(INDPRP('X2SPNORB')) = .TRUE.
98                     BQROP(INDPRP('Y2SPNORB')) = .TRUE.
99                     BQROP(INDPRP('Z2SPNORB')) = .TRUE.
100                  ELSE IF (LABEL .EQ. 'MNFSPNOR') THEN
101                     BQROP(INDPRP('X1MNF-SO')) = .TRUE.
102                     BQROP(INDPRP('Y1MNF-SO')) = .TRUE.
103                     BQROP(INDPRP('Z1MNF-SO')) = .TRUE.
104                  ELSE IF (LABEL .EQ. 'FERMI CO') THEN
105                     CALL GPOPEN(LUPROP,'AOPROPER','UNKNOWN',' ',
106     &                           'UNFORMATTED',IDUMMY,.FALSE.)
107 199                 CONTINUE
108                     IF (NXTLAB(LABFND,RTNLBL,LUPROP)) THEN
109                        IF (LABFND(1:3) .EQ. 'FC ')
110     &                       BQROP(INDPRP(LABFND)) = .TRUE.
111                        GOTO 199
112                     END IF
113                     CALL GPCLOSE(LUPROP,'KEEP')
114                  ELSE IF (LABEL .EQ. 'SPIN-DIP') THEN
115                     CALL GPOPEN(LUPROP,'AOPROPER','UNKNOWN',' ',
116     &                           'UNFORMATTED',IDUMMY,.FALSE.)
117 299                 CONTINUE
118                     IF (NXTLAB(LABFND,RTNLBL,LUPROP)) THEN
119                        IF (LABFND(1:3) .EQ. 'SD ')
120     &                       BQROP(INDPRP(LABFND)) = .TRUE.
121                        GOTO 299
122                     END IF
123                     CALL GPCLOSE(LUPROP,'KEEP')
124                  ELSE IF (LABEL(1:3) .EQ. 'PSO' .OR.
125     &                     LABEL(2:4) .EQ. 'PSO') THEN
126                     CALL GPOPEN(LUPROP,'AOPROPER','UNKNOWN',' ',
127     &                           'UNFORMATTED',IDUMMY,.FALSE.)
128 399                 CONTINUE
129                     IF (NXTLAB(LABFND,RTNLBL,LUPROP)) THEN
130                        IF (LABFND(1:3) .EQ. 'PSO')
131     &                       BQROP(INDPRP(LABFND)) = .TRUE.
132                        GOTO 399
133                     END IF
134                     CALL GPCLOSE(LUPROP,'KEEP')
135                  ELSE
136                     BQROP( INDPRP(LABEL)) = .TRUE.
137                     OPRANK( INDPRP(LABEL)) = IRANKB
138                  END IF
139               GO TO 100
140    2          CONTINUE
141                  BCFREQ = .TRUE.
142                  READ (LUCMD,*) NBQRFR
143                  IF (NBQRFR.LE.MBQRFR) THEN
144                     READ (LUCMD,*) (BQRFR(J),J=1,NBQRFR)
145                  ELSE
146                     WRITE (LUPRI,'(3(/,A,I5),/)')
147     *               ' NUMBER OF FREQUENCIES SPECIFIED    :',NBQRFR,
148     *               ' IS GREATER THAN THE ALLOWED NUMBER :',MBQRFR,
149     *               ' THE NUMBER IS RESET TO THE MAXIMUM :',MBQRFR
150                     READ (LUCMD,*) (BQRFR(J),J=1,MBQRFR),
151     *                              (FFFF,J=MBQRFR+1,NBQRFR)
152                     NBQRFR = MBQRFR
153                  END IF
154               GO TO 100
155    3          CONTINUE
156                  READ(LUCMD,'(BN,A,I8)')LABEL,IRANKC
157                  IF (LABEL .EQ. 'ANGMOM  ') THEN
158                     CQROP(INDPRP('XANGMOM ')) = .TRUE.
159                     CQROP(INDPRP('YANGMOM ')) = .TRUE.
160                     CQROP(INDPRP('ZANGMOM ')) = .TRUE.
161                  ELSE IF (LABEL .EQ. '1SPNORB ') THEN
162                     CQROP(INDPRP('X1SPNORB')) = .TRUE.
163                     CQROP(INDPRP('Y1SPNORB')) = .TRUE.
164                     CQROP(INDPRP('Z1SPNORB')) = .TRUE.
165                  ELSE IF (LABEL .EQ. '2SPNORB ') THEN
166                     CQROP(INDPRP('X2SPNORB')) = .TRUE.
167                     CQROP(INDPRP('Y2SPNORB')) = .TRUE.
168                     CQROP(INDPRP('Z2SPNORB')) = .TRUE.
169                  ELSE IF (LABEL .EQ. 'MNFSPNOR') THEN
170                     CQROP(INDPRP('X1MNF-SO')) = .TRUE.
171                     CQROP(INDPRP('Y1MNF-SO')) = .TRUE.
172                     CQROP(INDPRP('Z1MNF-SO')) = .TRUE.
173                  ELSE IF (LABEL .EQ. 'FERMI CO') THEN
174                     CALL GPOPEN(LUPROP,'AOPROPER','UNKNOWN',' ',
175     &                           'UNFORMATTED',IDUMMY,.FALSE.)
176 198                 CONTINUE
177                     IF (NXTLAB(LABFND,RTNLBL,LUPROP)) THEN
178                        IF (LABFND(1:3) .EQ. 'FC ')
179     &                       CQROP(INDPRP(LABFND)) = .TRUE.
180                        GOTO 198
181                     END IF
182                     CALL GPCLOSE(LUPROP,'KEEP')
183                  ELSE IF (LABEL .EQ. 'SPIN-DIP') THEN
184                     CALL GPOPEN(LUPROP,'AOPROPER','UNKNOWN',' ',
185     &                           'UNFORMATTED',IDUMMY,.FALSE.)
186 298                 CONTINUE
187                     IF (NXTLAB(LABFND,RTNLBL,LUPROP)) THEN
188                        IF (LABFND(1:3) .EQ. 'SD ')
189     &                       CQROP(INDPRP(LABFND)) = .TRUE.
190                        GOTO 298
191                     END IF
192                     CALL GPCLOSE(LUPROP,'KEEP')
193                  ELSE IF (LABEL(1:3) .EQ. 'PSO' .OR.
194     &                     LABEL(2:4) .EQ. 'PSO') THEN
195                     CALL GPOPEN(LUPROP,'AOPROPER','UNKNOWN',' ',
196     &                           'UNFORMATTED',IDUMMY,.FALSE.)
197 398                 CONTINUE
198                     IF (NXTLAB(LABFND,RTNLBL,LUPROP)) THEN
199                        IF (LABFND(1:3) .EQ. 'PSO')
200     &                       CQROP(INDPRP(LABFND)) = .TRUE.
201                        GOTO 398
202                     END IF
203                     CALL GPCLOSE(LUPROP,'KEEP')
204                  ELSE
205                     CQROP( INDPRP(LABEL)) = .TRUE.
206                     OPRANK( INDPRP(LABEL)) = IRANKC
207                  END IF
208               GO TO 100
209    4          CONTINUE
210                  BCFREQ = .TRUE.
211                  READ (LUCMD,*) NCQRFR
212                  IF (NCQRFR.LE.MCQRFR) THEN
213                     READ (LUCMD,*) (CQRFR(J),J=1,NCQRFR)
214                  ELSE
215                     WRITE (LUPRI,'(3(/,A,I5),/)')
216     *               ' NUMBER OF FREQUENCIES SPECIFIED    :',NCQRFR,
217     *               ' IS GREATER THAN THE ALLOWED NUMBER :',MCQRFR,
218     *               ' THE NUMBER IS RESET TO THE MAXIMUM :',MCQRFR
219                     READ (LUCMD,*) (CQRFR(J),J=1,MCQRFR),
220     *                              (FFFF,J=MCQRFR+1,NCQRFR)
221                     NCQRFR = MCQRFR
222                  END IF
223               GO TO 100
224    5          CONTINUE
225                  READ(LUCMD,'(BN,A,I8)')LABEL,IRANKA
226                  IF (LABEL .EQ. 'ANGMOM  ') THEN
227                     AQROP(INDPRP('XANGMOM ')) = .TRUE.
228                     AQROP(INDPRP('YANGMOM ')) = .TRUE.
229                     AQROP(INDPRP('ZANGMOM ')) = .TRUE.
230                  ELSE IF (LABEL .EQ. '1SPNORB ') THEN
231                     AQROP(INDPRP('X1SPNORB')) = .TRUE.
232                     AQROP(INDPRP('Y1SPNORB')) = .TRUE.
233                     AQROP(INDPRP('Z1SPNORB')) = .TRUE.
234                  ELSE IF (LABEL .EQ. '2SPNORB ') THEN
235                     AQROP(INDPRP('X2SPNORB')) = .TRUE.
236                     AQROP(INDPRP('Y2SPNORB')) = .TRUE.
237                     AQROP(INDPRP('Z2SPNORB')) = .TRUE.
238                  ELSE IF (LABEL .EQ. 'MNFSPNOR') THEN
239                     AQROP(INDPRP('X1MNF-SO')) = .TRUE.
240                     AQROP(INDPRP('Y1MNF-SO')) = .TRUE.
241                     AQROP(INDPRP('Z1MNF-SO')) = .TRUE.
242                  ELSE IF (LABEL .EQ. 'FERMI CO') THEN
243                     CALL GPOPEN(LUPROP,'AOPROPER','UNKNOWN',' ',
244     &                           'UNFORMATTED',IDUMMY,.FALSE.)
245 197                 CONTINUE
246                     IF (NXTLAB(LABFND,RTNLBL,LUPROP)) THEN
247                        IF (LABFND(1:3) .EQ. 'FC ')
248     &                       AQROP(INDPRP(LABFND)) = .TRUE.
249                        GOTO 197
250                     END IF
251                     CALL GPCLOSE(LUPROP,'KEEP')
252                  ELSE IF (LABEL .EQ. 'SPIN-DIP') THEN
253                     CALL GPOPEN(LUPROP,'AOPROPER','UNKNOWN',' ',
254     &                           'UNFORMATTED',IDUMMY,.FALSE.)
255 297                 CONTINUE
256                     IF (NXTLAB(LABFND,RTNLBL,LUPROP)) THEN
257                        IF (LABFND(1:3) .EQ. 'SD ')
258     &                       AQROP(INDPRP(LABFND)) = .TRUE.
259                        GOTO 297
260                     END IF
261                     CALL GPCLOSE(LUPROP,'KEEP')
262                  ELSE IF (LABEL(1:3) .EQ. 'PSO' .OR.
263     &                     LABEL(2:4) .EQ. 'PSO') THEN
264                     CALL GPOPEN(LUPROP,'AOPROPER','UNKNOWN',' ',
265     &                           'UNFORMATTED',IDUMMY,.FALSE.)
266 397                 CONTINUE
267                     IF (NXTLAB(LABFND,RTNLBL,LUPROP)) THEN
268                        IF (LABFND(1:3) .EQ. 'PSO')
269     &                       AQROP(INDPRP(LABFND)) = .TRUE.
270                        GOTO 397
271                     END IF
272                     CALL GPCLOSE(LUPROP,'KEEP')
273                  ELSE
274                     AQROP( INDPRP(LABEL)) = .TRUE.
275                     OPRANK(INDPRP(LABEL)) = IRANKA
276                  END IF
277               GO TO 100
278    6          CONTINUE
279                  READ (LUCMD,*) MAXITL
280               GO TO 100
281    7          CONTINUE
282                  READ (LUCMD,*) THCLR
283               GO TO 100
284    8          CONTINUE
285                  READ (LUCMD,*) MAXITO
286               GO TO 100
287    9          CONTINUE
288                  READ (LUCMD,*) IPRHYP
289               GO TO 100
290   10          CONTINUE
291                  LABEL='XDIPLEN'
292                  AQROP( INDPRP(LABEL)) = .TRUE.
293                  BQROP( INDPRP(LABEL)) = .TRUE.
294                  CQROP( INDPRP(LABEL)) = .TRUE.
295                  LABEL='YDIPLEN'
296                  AQROP( INDPRP(LABEL)) = .TRUE.
297                  BQROP( INDPRP(LABEL)) = .TRUE.
298                  CQROP( INDPRP(LABEL)) = .TRUE.
299                  LABEL='ZDIPLEN'
300                  AQROP( INDPRP(LABEL)) = .TRUE.
301                  BQROP( INDPRP(LABEL)) = .TRUE.
302                  CQROP( INDPRP(LABEL)) = .TRUE.
303               GO TO 100
304   11          CONTINUE
305                  E3TEST = .TRUE.
306               GO TO 100
307   12          CONTINUE
308               GO TO 100
309   13          CONTINUE
310                  READ (LUCMD,*) ISPINA,ISPINB,ISPINC
311               GO TO 100
312   14          CONTINUE
313                  A2TEST = .TRUE.
314               GO TO 100
315   15          CONTINUE
316                  X2TEST = .TRUE.
317               GO TO 100
318   16          CONTINUE
319                  REFCHK = .TRUE.
320               GO TO 100
321   17          CONTINUE
322                  READ(LUCMD,*) IAABB
323               GO TO 100
324   18          CONTINUE
325                  QRSPEC = .TRUE.
326                  QRSHG  = .TRUE.
327               GO TO 100
328   19          CONTINUE
329                  QRSPEC = .TRUE.
330                  QRPOCK = .TRUE.
331               GO TO 100
332   20          CONTINUE
333                  QRSPEC = .TRUE.
334                  QROPRF = .TRUE.
335               GO TO 100
336   21          CONTINUE
337                  COMFREQ = .TRUE.
338                  READ (LUCMD,*) NBQRFR
339                  IF (NBQRFR.LE.MBQRFR) THEN
340                     READ (LUCMD,*) (BQRFR(J),J=1,NBQRFR)
341                  ELSE
342                     WRITE (LUPRI,'(3(/,A,I5),/)')
343     *               ' NUMBER OF FREQUENCIES SPECIFIED    :',NBQRFR,
344     *               ' IS GREATER THAN THE ALLOWED NUMBER :',MBQRFR,
345     *               ' THE NUMBER IS RESET TO THE MAXIMUM :',MBQRFR
346                     READ (LUCMD,*) (BQRFR(J),J=1,MBQRFR),
347     *                              (FFFF,J=MBQRFR+1,NBQRFR)
348                     NBQRFR = MBQRFR
349                  END IF
350                  NCQRFR = NBQRFR
351                  DO J = 1,NBQRFR
352                     CQRFR(J) = BQRFR(J)
353                  END DO
354               GO TO 100
355   22          CONTINUE
356                  LABEL='XDIPLEN'
357                  AQROP( INDPRP(LABEL)) = .TRUE.
358                  BQROP( INDPRP(LABEL)) = .TRUE.
359                  CQROP( INDPRP(LABEL)) = .TRUE.
360               GO TO 100
361   23          CONTINUE
362                  LABEL='YDIPLEN'
363                  AQROP( INDPRP(LABEL)) = .TRUE.
364                  BQROP( INDPRP(LABEL)) = .TRUE.
365                  CQROP( INDPRP(LABEL)) = .TRUE.
366               GO TO 100
367   24          CONTINUE
368                  LABEL='ZDIPLEN'
369                  AQROP( INDPRP(LABEL)) = .TRUE.
370                  BQROP( INDPRP(LABEL)) = .TRUE.
371                  CQROP( INDPRP(LABEL)) = .TRUE.
372               GO TO 100
373   25          CONTINUE
374               READ(LUCMD,*)ISPINA
375               GO TO 100
376   26          CONTINUE
377               READ(LUCMD,*)ISPINB
378               GO TO 100
379   27          CONTINUE
380               READ(LUCMD,*)ISPINC
381               GO TO 100
382 28            CONTINUE
383                  SOCOLL = .TRUE.
384               GO TO 100
385 29            CONTINUE
386                  SSCOLL = .TRUE.
387               GO TO 100
388 30            CONTINUE
389                  QLOP = .TRUE.
390               GO TO 100
391            ELSE IF (PROMPT .EQ. '*') THEN
392               GO TO 300
393            ELSE
394               WRITE (LUPRI,'(/,3A,/)') ' PROMPT "',WORD,
395     *            '" NOT RECOGNIZED IN *QUADRA.'
396               CALL QUIT(' ILLEGAL PROMPT IN *QUADRA')
397            END IF
398         GO TO 100
399      END IF
400  300 CONTINUE
401      IF (THR_REDFAC .GT. 0.0D0) THEN
402         ICHANG = ICHANG + 1
403         WRITE (LUPRI,'(3A,1P,D10.2)') '@ INFO ',WORD1,
404     &   ' thresholds multiplied with general factor',THR_REDFAC
405         THCLR = THCLR*THR_REDFAC
406      END IF
407      NAQRTO = 0
408      NBQRTO = 0
409      NCQRTO = 0
410      IF (ICHANG .GT. 0) THEN
411         DO 500 I = 1,NPRLBL
412            IF (AQROP(I)) THEN
413                NAQRTO = NAQRTO + 1
414                ! For singlet states couple operator rank to excitation rank
415                IF (ISPIN.EQ.1) OPRANK(I) = ISPINA
416            END IF
417            IF (BQROP(I)) THEN
418                NBQRTO = NBQRTO + 1
419                IF (ISPIN.EQ.1) OPRANK(I) = ISPINB
420            END IF
421            IF (CQROP(I)) THEN
422                NCQRTO = NCQRTO + 1
423                IF (ISPIN.EQ.1) OPRANK(I) = ISPINC
424            END IF
425  500    CONTINUE
426C
427C  If special processes specified, we need to put the input freq. into
428C  the relevant frequency arrays.
429C
430         IF (QRSPEC) THEN
431            IF (BCFREQ) THEN
432               WRITE (LUPRI,'(/A/A/A)')
433     &  ' ERROR : inconsistent quadratic response input as',
434     &  ' .BFREQ and/or .CFREQ  has been specified together with'//
435     &     ' .SHG and/or .POCKEL',
436     &  ' Use only .FREQUE with .SHG, .POCKEL or .OPTREF!'
437               CALL QUIT('Inconsistent input for quadratic response')
438            END IF
439            NCQRFR = 0
440C
441            JFR = 0
442            IF (QRPOCK) THEN
443               NCQRFR = NCQRFR + 1
444               CQRFR(NCQRFR) = ZERO
445               DO IFR=1,NBQRFR
446                  IF (ABS(BQRFR(IFR)).LE.THRFRQ) JFR = IFR
447               ENDDO
448            ENDIF
449C
450            IF (QRSHG) THEN
451               DO IFR=1,NBQRFR
452                  IF (IFR.NE.JFR) THEN
453C                     IFR.EQ.JFR has already been included under QRPOCK
454                     NCQRFR = NCQRFR + 1
455                     CQRFR(NCQRFR) = BQRFR(IFR)
456                  END IF
457               END DO
458            ENDIF
459C
460            IF (QROPRF) THEN
461               DO IFR=1,NBQRFR
462                  IF (IFR.NE.JFR) THEN
463C                     IFR.EQ.JFR has already been included under QRPOCK
464                     NCQRFR = NCQRFR + 1
465                     CQRFR(NCQRFR) = -BQRFR(IFR)
466                  END IF
467               END DO
468            ENDIF
469C
470         ENDIF
471         IF (COMFREQ .AND. BCFREQ) THEN
472            WRITE (LUPRI,'(/A/A/A)')
473     &  ' *QUADRA ERROR : inconsistent quadratic response input as',
474     &  ' .BFREQ and/or .CFREQ  has been specified together with'//
475     &     ' .FREQUE',
476     &  ' Use either .FREQUE or .BFREQ/.CFREQ !!!!!'
477            CALL QUIT('Inconsistent input for quadratic response')
478         END IF
479      END IF
480      NQRCAL = MIN(NAQRTO,NBQRTO,NCQRTO,NBQRFR,NCQRFR)
481      IF  (HYPCAL .AND. NQRCAL.LE.0)  THEN
482        HYPCAL = .FALSE.
483        WRITE(LUPRI,'(/A)') ' Quadratic response input ignored because:'
484        IF (NAQRTO .EQ. 0) WRITE (LUPRI,'(A)')
485     *     ' - no A operators requested'
486        IF (NBQRTO .EQ. 0) WRITE (LUPRI,'(A)')
487     *     ' - no B operators requested'
488        IF (NBQRFR .LE. 0) WRITE (LUPRI,'(A,I8)')
489     &     ' - no B frequencies requested, NBQRFR =',NBQRFR
490        IF (NCQRTO .EQ. 0) WRITE (LUPRI,'(A)')
491     *     ' - no C operators requested'
492        IF (NCQRFR .LE. 0) WRITE (LUPRI,'(A,I8)')
493     &     ' - no C frequencies requested, NCQRFR =',NCQRFR
494      ELSE IF (HYPCAL) THEN
495         CALL HEADER('Quadratic Response calculation',0)
496         WRITE (LUPRI,'(A,L2)')
497     *      ' First hyperpolarizability calculation : HYPCAL='
498     *      ,HYPCAL
499         WRITE (LUPRI,'(/A,I5)') ' Spin of operator A , ISPINA=',ISPINA
500         WRITE (LUPRI,'( A,I5)') ' Spin of operator B , ISPINB=',ISPINB
501         WRITE (LUPRI,'( A,I5)') ' Spin of operator C , ISPINC=',ISPINC
502         IF (ISPINA+ISPINB+ISPINC .GT. 0) THEN
503            TRPLET = .TRUE.
504            TRPFLG = .TRUE.
505         END IF
506         WRITE(LUPRI,'(/I3,A,(1P,5D14.6))')
507     *      NBQRFR,' B-frequencies', (BQRFR(I),I=1,NBQRFR)
508         WRITE(LUPRI,'(I3,A,(1P,5D14.6))')
509     *      NCQRFR,' C-frequencies', (CQRFR(I),I=1,NCQRFR)
510         IF (SOCOLL) WRITE (LUPRI,'(/A)')
511     *        ' Quadratic response functions will be collected as '//
512     *        'SO-corrections to nuclear shieldings'
513         IF (SSCOLL) WRITE (LUPRI,'(/A)')
514     *        ' Quadratic response functions will be collected as '//
515     *        'SO-corrections to reduced spin-spin couplings'
516         IF (FLAG(16) .AND. .NOT.INERSI) THEN
517            WRITE(LUPRI,'(/A,L2)')
518     &      ' Equilibrium solvent model requested            : SOLVNT ='
519     &      ,FLAG(16)
520            WRITE(LUPRI,'(A,F8.4)')
521     &      '    Dielectric constant                         : EPSOL  ='
522     &      ,EPSOL
523         END IF
524         IF (FLAG(16) .AND. INERSI) THEN
525            WRITE(LUPRI,'(/A,L2)')
526     &      ' Non-equilibrium solvent model requested        : INERSI ='
527     &      ,INERSI
528            WRITE(LUPRI,'(A,F8.4)')
529     &      '    Static dielectric constant                  : EPSTAT ='
530     &      ,EPSTAT
531            WRITE(LUPRI,'(A,F8.4)')
532     &      '    Optical dielectric constant                 : EPSOL  ='
533     &      ,EPSOL
534         END IF
535
536         IF (QM3) WRITE(LUPRI,'(/A)')
537     &      ' QM/MM response calculation (.QM3)'
538
539         IF (QMMM) WRITE(LUPRI,'(/A)')
540     &      ' QM/MM response calculation (.QMMM)'
541
542         IF (USE_PELIB()) THEN
543            WRITE(LUPRI,'(/A)')
544     &         'Environment effects are included (.PELIB)'
545         END IF
546
547         IF (QMNPMM) THEN
548            WRITE(LUPRI,'(/A)') ' Heterogeneous environment is '//
549     &                          ' modeled using the QM/NP/MM scheme'
550         END IF
551
552         WRITE(LUPRI,'()')
553
554         IF (QRSHG) WRITE (LUPRI,'(A)')
555     &      ' Second harmonic generation;        SHG'
556         IF (QRPOCK) WRITE (LUPRI,'(A)')
557     &      ' Electro-optical Pockels effect;    POCKEL'
558         IF (QROPRF) WRITE (LUPRI,'(A)')
559     &      ' Optical refractivity calculation;  OPTREF'
560
561         WRITE(LUPRI,'(/A,I4)')
562     *      ' Print level                                    : IPRHYP ='
563     *       ,IPRHYP
564         WRITE(LUPRI,'(A,I4)')
565     *      ' Maximum number of iterations in lin.rsp. solver: MAXITL ='
566     *       ,MAXITL
567         WRITE(LUPRI,'(A,1P,D10.3)')
568     *      ' Threshold for convergence of linear resp. eq.s : THCLR  ='
569     *      ,THCLR
570         WRITE(LUPRI,'(A,I4)')
571     *      ' Maximum iterations in optimal orbital algorithm: MAXITO ='
572     *      ,MAXITO
573         IF (DIROIT) WRITE (LUPRI,'(A,L2)')
574     *      ' Direct one-index transformation                : DIROIT ='
575     *      ,DIROIT
576         IF (E3TEST) WRITE (LUPRI,'(/A)')
577     *      ' Test of contributions from E3 and S3 terms.'
578         IF (A2TEST) WRITE (LUPRI,'(/A)')
579     *      ' Test of contributions from A2 terms.'
580         IF (X2TEST) WRITE (LUPRI,'(/A)')
581     *      ' Test of contributions from X2 terms.'
582      END IF
583C
584C *** END OF HYPINP
585C
586      RETURN
587      END
588C  /* Deck smoinp */
589      SUBROUTINE SMOINP(WORD)
590C
591C Sonia: added MCDBTE section (1999)
592C
593      use pelib_interface, only: use_pelib
594#include "implicit.h"
595C
596#include "priunit.h"
597#include "infrsp.h"
598#include "inforb.h"
599#include "rspprp.h"
600#include "infsmo.h"
601#include "inflr.h"
602#include "infpp.h"
603#include "infpri.h"
604#include "infspi.h"
605#include "inflin.h"
606#include "infhso.h"
607#include "maxorb.h"
608#include "infinp.h"
609#include "inftap.h"
610! gnrinf.h: QM3, THR_REDFAC
611#include "gnrinf.h"
612C
613      LOGICAL NEWDEF
614      PARAMETER ( NTABLE = 33)                                          !MK
615      CHARACTER PROMPT*1, WORD*7, TABLE(NTABLE)*7, WORD1*7
616      CHARACTER*8 LABEL
617C
618      DATA TABLE /'.BPROP ', '.BFREQ ', '.APROP ', '.ROOTS ',
619     *            '.MAXITP', '.THCLR ', '.MAXITL', '.THCPP ',
620     *            '.MAXITO', '.PRINT ', '.DIPLEN', '.E3TEST',
621     *            '.ISPABC', '.TPCD  ', '.X2TEST', '.A2TEST',
622     *            '.DIPLNX', '.DIPLNY', '.DIPLNZ', '.FREQUE',
623     *            '.SINGLE', '.PHOSPH', '.MNFPHO', '.TWO-PH',
624     *            '.MCDBTE', '.MCDPRJ', '.ECPHOS', '.PHOSPV',           !MK
625     *            '.DIPVEL', '.CPPHOL', '.CPPHOV', '.CPPHMF',           !MK
626     *            '.CPPHEC'/                                            !MK
627C
628C READ IN  INPUT
629C
630      NEWDEF = (WORD .EQ. '*QUADRA')
631      ICHANG = 0
632      IF (NEWDEF) THEN
633         SOMOM = .TRUE.
634         JSPABC = 0
635         WORD1 = WORD
636  100    CONTINUE
637            READ (LUCMD, '(A7)') WORD
638            CALL UPCASE(WORD)
639            PROMPT = WORD(1:1)
640            IF (PROMPT .EQ. '!' .OR. PROMPT .EQ. '#') GO TO 100
641            IF (PROMPT .EQ. '.') THEN
642               ICHANG = ICHANG + 1
643               DO 200 I = 1, NTABLE
644                  IF (TABLE(I) .EQ. WORD) THEN
645                     GO TO (1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16,
646     &                 17,18,19,20,21,22,23,24,25,26,27,28,29,30,31,
647     &                 32,33), I                                        !MK
648                  END IF
649  200          CONTINUE
650               IF (WORD .EQ. '.OPTION') THEN
651                 CALL PRTAB(NTABLE,TABLE,WORD1//' input keywords',LUPRI)
652                 GO TO 100
653               END IF
654               WRITE (LUPRI,'(/,3A,/)') ' KEYWORD "',WORD,
655     *            '" NOT RECOGNIZED IN SMOINP.'
656               CALL PRTAB(NTABLE,TABLE,WORD1//' input keywords',LUPRI)
657               CALL QUIT(' ILLEGAL KEYWORD IN SMOINP ')
658 1             CONTINUE
659                  READ(LUCMD,'( A )')LABEL
660                  BSMOP( INDPRP(LABEL)) = .TRUE.
661               GO TO 100
662 2             CONTINUE
663                  READ (LUCMD,*) NBSMFR
664                  IF (NBSMFR.LE.MBSMFR) THEN
665                     READ (LUCMD,*) (BSMFR(J),J=1,NBSMFR)
666                  ELSE
667                     WRITE (LUPRI,'(3(/,A,I5),/)')
668     *               ' NUMBER OF B FREQUENCIES SPECIFIED  :',NBSMFR,
669     *               ' IS GREATER THAN THE ALLOWED NUMBER :',MBSMFR,
670     *               ' THE NUMBER IS RESET TO THE MAXIMUM :',MBSMFR
671                     READ (LUCMD,*) (BSMFR(J),J=1,MBSMFR),
672     *                              (FFFF,J=MBSMFR+1,NBSMFR)
673                     NBSMFR = MBSMFR
674                  END IF
675               GO TO 100
676 3             CONTINUE
677                  READ(LUCMD,'( A )')LABEL
678                  ASMOP( INDPRP(LABEL)) = .TRUE.
679               GO TO 100
680 4             CONTINUE
681                  READ (LUCMD,*) (NSMCNV(MULD2H(J,LSYMRF)),J=1,NSYM)
682               GO TO 100
683 5             CONTINUE
684                  READ (LUCMD,*) MAXITL
685               GO TO 100
686 6             CONTINUE
687                  READ (LUCMD,*) THCLR
688               GO TO 100
689 7             CONTINUE
690                  READ (LUCMD,*) MAXITP
691               GO TO 100
692 8             CONTINUE
693                  READ (LUCMD,*) THCPP
694               GO TO 100
695 9             CONTINUE
696                  READ (LUCMD,*) MAXITO
697               GO TO 100
698 10            CONTINUE
699                  READ (LUCMD,*) IPRSMO
700               GO TO 100
701 11            CONTINUE
702                  LABEL='XDIPLEN'
703                  ASMOP( INDPRP(LABEL)) = .TRUE.
704                  BSMOP( INDPRP(LABEL)) = .TRUE.
705                  LABEL='YDIPLEN'
706                  ASMOP( INDPRP(LABEL)) = .TRUE.
707                  BSMOP( INDPRP(LABEL)) = .TRUE.
708                  LABEL='ZDIPLEN'
709                  ASMOP( INDPRP(LABEL)) = .TRUE.
710                  BSMOP( INDPRP(LABEL)) = .TRUE.
711               GO TO 100
712 29            CONTINUE                                                 !MK
713                  LABEL='XDIPVEL'                                       !MK
714                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
715                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
716                  LABEL='YDIPVEL'                                       !MK
717                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
718                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
719                  LABEL='ZDIPVEL'                                       !MK
720                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
721                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
722               GO TO 100                                                !MK
723 12            CONTINUE
724                  E3TEST = .TRUE.
725               GO TO 100
726 13            CONTINUE
727                  JSPABC = JSPABC + 1
728                  READ (LUCMD,*) ISPINA,ISPINB,ISPINC
729               GO TO 100
730C              Option 14 .TPACD calculates tensor components for Tinoco
731C              original, gauge invariant formulation
732 14            CONTINUE
733                  TWOPHO = .TRUE.
734                  LTPCD  = .TRUE.
735C                 A operators
736                  LABEL='XDIPVEL'
737                  ASMOP( INDPRP(LABEL)) = .TRUE.
738                  LABEL='YDIPVEL'
739                  ASMOP( INDPRP(LABEL)) = .TRUE.
740                  LABEL='ZDIPVEL'
741                  ASMOP( INDPRP(LABEL)) = .TRUE.
742
743                  LABEL='XANGMOM'
744                  ASMOP( INDPRP(LABEL)) = .TRUE.
745                  LABEL='YANGMOM'
746                  ASMOP( INDPRP(LABEL)) = .TRUE.
747                  LABEL='ZANGMOM'
748                  ASMOP( INDPRP(LABEL)) = .TRUE.
749
750                  LABEL='XXROTSTR'
751                  ASMOP( INDPRP(LABEL)) = .TRUE.
752                  LABEL='YYROTSTR'
753                  ASMOP( INDPRP(LABEL)) = .TRUE.
754                  LABEL='ZZROTSTR'
755                  ASMOP( INDPRP(LABEL)) = .TRUE.
756                  LABEL='XYROTSTR'
757                  ASMOP( INDPRP(LABEL)) = .TRUE.
758                  LABEL='XZROTSTR'
759                  ASMOP( INDPRP(LABEL)) = .TRUE.
760                  LABEL='YZROTSTR'
761                  ASMOP( INDPRP(LABEL)) = .TRUE.
762
763C                 B operators
764
765                  LABEL='XDIPVEL'
766                  BSMOP( INDPRP(LABEL)) = .TRUE.
767                  LABEL='YDIPVEL'
768                  BSMOP( INDPRP(LABEL)) = .TRUE.
769                  LABEL='ZDIPVEL'
770                  BSMOP( INDPRP(LABEL)) = .TRUE.
771               GO TO 100
772 15            CONTINUE
773                  X2TEST = .TRUE.
774               GO TO 100
775 16            CONTINUE
776                  A2TEST = .TRUE.
777               GO TO 100
778 17            CONTINUE
779                  LABEL='XDIPLEN'
780                  ASMOP( INDPRP(LABEL)) = .TRUE.
781                  BSMOP( INDPRP(LABEL)) = .TRUE.
782               GO TO 100
783 18            CONTINUE
784                  LABEL='YDIPLEN'
785                  ASMOP( INDPRP(LABEL)) = .TRUE.
786                  BSMOP( INDPRP(LABEL)) = .TRUE.
787               GO TO 100
788 19            CONTINUE
789                  LABEL='ZDIPLEN'
790                  ASMOP( INDPRP(LABEL)) = .TRUE.
791                  BSMOP( INDPRP(LABEL)) = .TRUE.
792               GO TO 100
793 20            CONTINUE
794               GO TO 2
795 21            CONTINUE
796               GO TO 100
797 22            CONTINUE
798                  PHOSPH = .TRUE.
799                  LABEL='XDIPLEN '
800                  ASMOP( INDPRP(LABEL)) = .TRUE.
801                  LABEL='YDIPLEN '
802                  ASMOP( INDPRP(LABEL)) = .TRUE.
803                  LABEL='ZDIPLEN '
804                  ASMOP( INDPRP(LABEL)) = .TRUE.
805                  LABEL='X SPNORB'
806                  BSMOP( INDPRP(LABEL)) = .TRUE.
807                  LABEL='Y SPNORB'
808                  BSMOP( INDPRP(LABEL)) = .TRUE.
809                  LABEL='Z SPNORB'
810                  BSMOP( INDPRP(LABEL)) = .TRUE.
811               GO TO 100
812 28            CONTINUE                                                 !MK
813                  PHOSPV = .TRUE.                                       !MK
814                  LABEL='XDIPVEL '                                      !MK
815                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
816                  LABEL='YDIPVEL '                                      !MK
817                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
818                  LABEL='ZDIPVEL '                                      !MK
819                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
820                  LABEL='X SPNORB'                                      !MK
821                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
822                  LABEL='Y SPNORB'                                      !MK
823                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
824                  LABEL='Z SPNORB'                                      !MK
825                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
826               GO TO 100                                                !MK
827 30            CONTINUE                                                 !MK
828                  CPPHOL = .TRUE.                                       !MK
829                  LABEL='XDIPLEN '                                      !MK
830                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
831                  LABEL='YDIPLEN '                                      !MK
832                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
833                  LABEL='ZDIPLEN '                                      !MK
834                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
835                  LABEL='X SPNORB'                                      !MK
836                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
837                  LABEL='Y SPNORB'                                      !MK
838                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
839                  LABEL='Z SPNORB'                                      !MK
840                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
841                  LABEL='XANGMOM'                                       !MK
842                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
843                  LABEL='YANGMOM'                                       !MK
844                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
845                  LABEL='ZANGMOM'                                       !MK
846                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
847               GO TO 100                                                !MK
848 31            CONTINUE                                                 !MK
849                  CPPHOV = .TRUE.                                       !MK
850                  LABEL='XDIPVEL '                                      !MK
851                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
852                  LABEL='YDIPVEL '                                      !MK
853                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
854                  LABEL='ZDIPVEL '                                      !MK
855                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
856                  LABEL='X SPNORB'                                      !MK
857                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
858                  LABEL='Y SPNORB'                                      !MK
859                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
860                  LABEL='Z SPNORB'                                      !MK
861                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
862                  LABEL='XANGMOM'                                       !MK
863                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
864                  LABEL='YANGMOM'                                       !MK
865                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
866                  LABEL='ZANGMOM'                                       !MK
867                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
868               GO TO 100                                                !MK
869 32            CONTINUE                                                 !MK
870                  CPPHMF = .TRUE.                                       !MK
871                  LABEL='XDIPVEL '                                      !MK
872                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
873                  LABEL='YDIPVEL '                                      !MK
874                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
875                  LABEL='ZDIPVEL '                                      !MK
876                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
877                  LABEL='XANGMOM'                                       !MK
878                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
879                  LABEL='YANGMOM'                                       !MK
880                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
881                  LABEL='ZANGMOM'                                       !MK
882                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
883                  LABEL='X1MNF-SO'                                      !MK
884                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
885                  LABEL='Y1MNF-SO'                                      !MK
886                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
887                  LABEL='Z1MNF-SO'                                      !MK
888                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
889               GO TO 100                                                !MK
890 33            CONTINUE                                                 !MK
891                  CPPHEC = .TRUE.                                       !MK
892                  LABEL='XDIPVEL '                                      !MK
893                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
894                  LABEL='YDIPVEL '                                      !MK
895                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
896                  LABEL='ZDIPVEL '                                      !MK
897                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
898                  LABEL='X1SPNSCA'                                      !MK
899                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
900                  LABEL='Y1SPNSCA'                                      !MK
901                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
902                  LABEL='Z1SPNSCA'                                      !MK
903                  BSMOP( INDPRP(LABEL)) = .TRUE.                        !MK
904                  LABEL='XANGMOM'                                       !MK
905                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
906                  LABEL='YANGMOM'                                       !MK
907                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
908                  LABEL='ZANGMOM'                                       !MK
909                  ASMOP( INDPRP(LABEL)) = .TRUE.                        !MK
910               GO TO 100                                                !MK
911 23            CONTINUE
912                  MNFPHO = .TRUE.
913                  LABEL='XDIPLEN '
914                  ASMOP( INDPRP(LABEL)) = .TRUE.
915                  LABEL='YDIPLEN '
916                  ASMOP( INDPRP(LABEL)) = .TRUE.
917                  LABEL='ZDIPLEN '
918                  ASMOP( INDPRP(LABEL)) = .TRUE.
919                  LABEL='X1MNF-SO'
920                  BSMOP( INDPRP(LABEL)) = .TRUE.
921                  LABEL='Y1MNF-SO'
922                  BSMOP( INDPRP(LABEL)) = .TRUE.
923                  LABEL='Z1MNF-SO'
924                  BSMOP( INDPRP(LABEL)) = .TRUE.
925               GO TO 100
926 27            CONTINUE
927                  ECPHOS = .TRUE.
928                  LABEL='XDIPLEN '
929                  ASMOP( INDPRP(LABEL)) = .TRUE.
930                  LABEL='YDIPLEN '
931                  ASMOP( INDPRP(LABEL)) = .TRUE.
932                  LABEL='ZDIPLEN '
933                  ASMOP( INDPRP(LABEL)) = .TRUE.
934                  LABEL='X1SPNSCA'
935                  BSMOP( INDPRP(LABEL)) = .TRUE.
936                  LABEL='Y1SPNSCA'
937                  BSMOP( INDPRP(LABEL)) = .TRUE.
938                  LABEL='Z1SPNSCA'
939                  BSMOP( INDPRP(LABEL)) = .TRUE.
940               GO TO 100
941 24            CONTINUE
942                  TWOPHO = .TRUE.
943                  LABEL='XDIPLEN'
944                  ASMOP( INDPRP(LABEL)) = .TRUE.
945                  BSMOP( INDPRP(LABEL)) = .TRUE.
946                  LABEL='YDIPLEN'
947                  ASMOP( INDPRP(LABEL)) = .TRUE.
948                  BSMOP( INDPRP(LABEL)) = .TRUE.
949                  LABEL='ZDIPLEN'
950                  ASMOP( INDPRP(LABEL)) = .TRUE.
951                  BSMOP( INDPRP(LABEL)) = .TRUE.
952               GO TO 100
953 25            CONTINUE
954                  MCDCAL=.TRUE.
955                  LABEL='XDIPLEN'
956                  ASMOP( INDPRP(LABEL)) = .TRUE.
957                  LABEL='YDIPLEN'
958                  ASMOP( INDPRP(LABEL)) = .TRUE.
959                  LABEL='ZDIPLEN'
960                  ASMOP( INDPRP(LABEL)) = .TRUE.
961                  LABEL='XANGMOM'
962                  BSMOP( INDPRP(LABEL)) = .TRUE.
963                  LABEL='YANGMOM'
964                  BSMOP( INDPRP(LABEL)) = .TRUE.
965                  LABEL='ZANGMOM'
966                  BSMOP( INDPRP(LABEL)) = .TRUE.
967               GO TO 100
968 26            CONTINUE
969                  !sonia: mcd B term calculation
970                  !with projection of the redundant exc state vector
971                  MCDPRJ=.TRUE.
972                  MCDCAL=.TRUE.
973                  LABEL='XDIPLEN'
974                  ASMOP( INDPRP(LABEL)) = .TRUE.
975                  LABEL='YDIPLEN'
976                  ASMOP( INDPRP(LABEL)) = .TRUE.
977                  LABEL='ZDIPLEN'
978                  ASMOP( INDPRP(LABEL)) = .TRUE.
979                  LABEL='XANGMOM'
980                  BSMOP( INDPRP(LABEL)) = .TRUE.
981                  LABEL='YANGMOM'
982                  BSMOP( INDPRP(LABEL)) = .TRUE.
983                  LABEL='ZANGMOM'
984                  BSMOP( INDPRP(LABEL)) = .TRUE.
985               GO TO 100
986            ELSE IF (PROMPT .EQ. '*') THEN
987               GO TO 300
988            ELSE
989               WRITE (LUPRI,'(/,3A,/)') ' PROMPT "',WORD,
990     *            '" NOT RECOGNIZED IN SMOINP.'
991               CALL QUIT(' ILLEGAL PROMPT IN SMOINP ')
992            END IF
993         GO TO 100
994      END IF
995  300 CONTINUE
996      IF (THR_REDFAC .GT. 0.0D0) THEN
997         ICHANG = ICHANG + 1
998         WRITE (LUPRI,'(3A,1P,D10.2)') '@ INFO ',WORD1,
999     &   ' thresholds multiplied with general factor',THR_REDFAC
1000         THCLR = THCLR*THR_REDFAC
1001         THCPP = THCPP*THR_REDFAC
1002      END IF
1003C
1004C If two-photon calculation then we assume that the frequencies of the
1005C photons match the excitation. At this point we set the frequencies
1006C to dummy values, later they become half of the excitation energies.
1007C
1008      IF (TWOPHO) THEN
1009         NSMCNV_TOT = ISUM(NSYM,NSMCNV,1)
1010         IF (NSMCNV_TOT .EQ. 0) THEN
1011            WRITE (LUPRI,'(/2A)')
1012     * ' Two-photon input ignored because no excited states specified.',
1013     * ' Use the .ROOTS keyword.'
1014            CALL QUIT('***** Input error for two-photon *****')
1015         END IF
1016         NBSMFR = NSMCNV_TOT
1017         DO J=1,NBSMFR
1018            BSMFR(J)=-1.0D0
1019         END DO
1020      END IF
1021C
1022      NASMTO = 0
1023      NBSMTO = 0
1024      IF (ICHANG .GT. 0) THEN
1025         DO 500 I = 1,NPRLBL
1026            IF (ASMOP(I)) NASMTO = NASMTO + 1
1027            IF (BSMOP(I)) NBSMTO = NBSMTO + 1
1028 500     CONTINUE
1029      END IF
1030      NSMCAL = MIN(NASMTO,NBSMTO,NBSMFR)
1031      IF  (SOMOM .AND. NSMCAL.LE.0)  THEN
1032         SOMOM = .FALSE.
1033         WRITE (LUPRI,'(/A)')
1034     & ' Quadratic Response single residue calculation ignored because:'
1035         IF (NASMTO .EQ. 0) WRITE (LUPRI,'(A)')
1036     &      ' - no A operators requested'
1037         IF (NBSMTO .EQ. 0) WRITE (LUPRI,'(A)')
1038     &      ' - no B operators requested'
1039         IF (NBSMFR .LE. 0) WRITE (LUPRI,'(A,I8)')
1040     &      ' - no B frequencies requested, NBSMFR =',NBSMFR
1041      ELSE IF (SOMOM) THEN
1042         CALL HEADER('Quadratic Response single residue calculation',0)
1043         WRITE(LUPRI,'(A)')
1044     &   ' That is, calculation of second order moments.'
1045         IF (PHOSPH) THEN
1046           WRITE(LUPRI,'(/A)') ' Phosphorescence calculation requested'
1047     &     //' with full spin-orbit operator.'
1048         END IF
1049         IF (MNFPHO) THEN
1050           WRITE(LUPRI,'(/A)') ' Phosphorescence calculation requested'
1051     &     //' with atomic mean-field spin-orbit operator (AMFI).'
1052         END IF
1053         IF (ECPHOS) THEN
1054           WRITE(LUPRI,'(/A)') ' Phosphorescence calculation requested'
1055     &     //' with effective charge spin-orbit operator.'
1056         END IF
1057         IF (PHOSPH .OR. MNFPHO .OR. ECPHOS .OR. PHOSPV .OR.            !MK
1058     &       CPPHOL .OR. CPPHOV .OR. CPPHMF .OR. CPPHEC) THEN           !MK
1059                 IF (JSPABC .GT. 0 .AND. ISPINA .NE. 0 .AND.
1060     &          ISPINB .NE. 1 .AND. ISPINC .NE. 1) THEN
1061               NWARN = NWARN + 1
1062               WRITE (LUPRI,'(/A/A,3I4/A)')
1063     &         ' WARNING : inconsistent phosphorescence input as',
1064     &         ' .ISPABC has been specified with',ISPINA,ISPINB,ISPINC,
1065     &         ' .ISPABC is reset to               0   1   1'
1066            END IF
1067            ISPINA = 0
1068            ISPINB = 1
1069            ISPINC = 1
1070C
1071            I = 0
1072            DO J = 1,NSYM
1073            IF (NSMCNV(J) .GT. MXPHOS) THEN
1074               I = I + 1
1075               NSMCNV(J) = MXPHOS
1076            END IF
1077            END DO
1078            IF (I .GT. 0) THEN
1079               NWARN = NWARN + 1
1080               WRITE (LUPRI,'(/A,I3,A/I2,A/)')
1081     &         ' WARNING: infhso.h is currently dimensioned to max.',
1082     &            MXPHOS,' excited states of any symmetry.',
1083     &         I,' specification(s) under .ROOTS have been'//
1084     &            ' reduced to this value.'
1085            END IF
1086         END IF
1087         IF (TWOPHO) THEN
1088            WRITE(LUPRI,'(/A/A)')
1089     &      ' Two-photon transition processes computed (.TWO-PHOTON).',
1090     &      ' Both photon frequencies are half the excitation energy.'
1091         ELSE
1092            WRITE(LUPRI,'(/I3,A,(1P,5D14.6))')
1093     &      NBSMFR,' B-frequencies', (BSMFR(I),I=1,NBSMFR)
1094         END IF
1095         WRITE (LUPRI,'(/A,I5)') ' Spin of operator A , ISPINA=',ISPINA
1096         WRITE (LUPRI,'( A,I5)') ' Spin of operator B , ISPINB=',ISPINB
1097         WRITE (LUPRI,'(2A,I5)') ' Spin of operator C ,'
1098     &    ,' (Excitation energy) ISPINC=',ISPINC
1099         IF (ISPINA+ISPINB+ISPINC .GT. 0) THEN
1100            TRPLET = .TRUE.
1101            TRPFLG = .TRUE.
1102         END IF
1103         IF (FLAG(16) .AND. .NOT.INERSI) THEN
1104            WRITE(LUPRI,'(/A,L2)')
1105     &      ' Equilibrium solvent model requested            : SOLVNT ='
1106     &      ,FLAG(16)
1107            WRITE(LUPRI,'(A,F8.4)')
1108     &      '    Dielectric constant                         : EPSOL  ='
1109     &      ,EPSOL
1110         END IF
1111         IF (FLAG(16) .AND. INERSI) THEN
1112            WRITE(LUPRI,'(/A,L2)')
1113     &      ' Non-equilibrium solvent model requested        : INERSI ='
1114     &      ,INERSI
1115            WRITE(LUPRI,'(A,F8.4)')
1116     &      '    Static dielectric constant                  : EPSTAT ='
1117     &      ,EPSTAT
1118            WRITE(LUPRI,'(A,F8.4)')
1119     &      '    Optical dielectric constant                 : EPSOL  ='
1120     &      ,EPSOL
1121         END IF
1122
1123         IF (QM3) WRITE(LUPRI,'(/A)')
1124     &      ' QM/MM response calculation (.QM3)'
1125
1126         IF (QMMM) WRITE(LUPRI,'(/A)')
1127     &      ' QM/MM response calculation (.QMMM)'
1128
1129         IF (USE_PELIB()) THEN
1130            WRITE(LUPRI,'(/A)')
1131     &         'Environment effects are included (.PELIB)'
1132         END IF
1133
1134         IF (MCDCAL) THEN
1135            WRITE (LUPRI,'(A,L1)')
1136     *       ' B of Magnetic Circular Dichroism requested    : MCDCAL ='
1137     *      ,MCDCAL
1138            !sonia new
1139            if (MCDPRJ) then
1140               WRITE (LUPRI,'(A,L1)')
1141     *       ' Projected B of Magnetic Circular Dichroism    : MCDPRJ='
1142     *        ,MCDPRJ
1143            end if
1144         END IF
1145         WRITE(LUPRI,'(/A,I5)')
1146     *      ' Print level                                    : IPRSMO ='
1147     *      ,IPRSMO
1148         WRITE(LUPRI,'(A,1P,D10.3)')
1149     *      ' Threshold for convergence in RSPPP             : THCPP  ='
1150     *      ,THCPP
1151         WRITE(LUPRI,'(A,I5)')
1152     *      ' Maximum number of iterations in RSPPP          : MAXITP ='
1153     *      ,MAXITP
1154         WRITE(LUPRI,'(A,1P,D10.3)')
1155     *      ' Threshold for convergence in RSPLR             : THCLR  ='
1156     *      ,THCLR
1157         WRITE(LUPRI,'(A,I5)')
1158     *      ' Maximum number of iterations in RSPLR          : MAXITL ='
1159     *      ,MAXITL
1160         WRITE(LUPRI,'(A,I5)')
1161     *      ' Maximum iterations in optimal orbital algorithm: MAXITO ='
1162     *      ,MAXITO
1163         IF (E3TEST) WRITE (LUPRI,'(/A)')
1164     *      ' Test of contributions from E3 and S3 terms.'
1165      END IF
1166C
1167C *** END OF SMOINP
1168C
1169      RETURN
1170      END
1171C  /* Deck exiind */
1172      SUBROUTINE EXIIND(NEVAL)
1173C
1174C SET UP A LIST WITH THE NUMBER OF EIGENVECTORS THAT HAS
1175C TO BE SOLVED
1176C
1177#include "implicit.h"
1178C
1179#include "rspprp.h"
1180#include "indqr.h"
1181#include "inforb.h"
1182#include "infpri.h"
1183#include "infsmo.h"
1184C
1185      DIMENSION NEVAL(*)
1186C
1187C If two-photon calculation is specified we only compute a restricted
1188C number of excitation vectors.
1189C
1190      IF (TWOPHO) THEN
1191         DO ICSYM = 1,NSYM
1192         DO IBSYM = 1,NSYM
1193            IASYM = MULD2H(ICSYM,IBSYM)
1194            IF ( (NSMCNV(ICSYM) .GT. 0) .AND. (NBSMOP(IBSYM) .GT. 0)
1195     *          .AND. (NASMOP(IASYM).GT.0) ) THEN
1196               DO IC = 1,NSMCNV(ICSYM)
1197                  INUM = INQREX(ICSYM,IC)
1198               END DO
1199            END IF
1200         END DO
1201         END DO
1202      ELSE
1203C
1204C The general case computes all combinations.
1205C
1206         DO 100 ISYM = 1,NSYM
1207            DO 200 IROOT = 1,NEVAL(ISYM)
1208               INUM = INQREX(ISYM,IROOT)
1209 200        CONTINUE
1210 100     CONTINUE
1211C
1212      END IF
1213C
1214      RETURN
1215      END
1216C  /* Deck inqrex */
1217      INTEGER FUNCTION INQREX(ISYM,IROOT)
1218C
1219C FIND THE NUMBER THAT HAS BEEN GIVEN TO THE INUM EXCITATION OF
1220C SYMMETRY ISYM
1221C
1222#include "implicit.h"
1223C
1224#include "priunit.h"
1225#include "infrsp.h"
1226#include "rspprp.h"
1227#include "indqr.h"
1228#include "infpri.h"
1229C
1230      DO 100 I=1,NEXLBL
1231         IF ( (ISYM .EQ. ISEXQR(I) ) .AND. ( IROOT.EQ.JEXQR(I))) THEN
1232            INQREX = I
1233            RETURN
1234         END IF
1235 100  CONTINUE
1236      NEXLBL = NEXLBL + 1
1237      IF ( NEXLBL.GT.MXEXQR ) THEN
1238         WRITE(LUPRI,'(/A,I6)')
1239     &      'ERROR INQREX: Number of QR excitation operators exceeds'
1240     &      //' the allowed number, MXEXQR =',MXEXQR
1241         CALL QUIT(' INQREX: TOO MANY EXCITATION OPERATORS')
1242      END IF
1243      JEXQR(NEXLBL) = IROOT
1244      ISEXQR(NEXLBL) = ISYM
1245      INQREX = NEXLBL
1246      IF (IPRRSP.GT.10) THEN
1247         WRITE(LUPRI,'(/A)')' INQREX: NEW EXCITATION '
1248         WRITE(LUPRI,'(A,I5,2X,A,I5,2X,A,I5)')
1249     *   ' EXCITATION NUMBER:',NEXLBL,'SYMMETRY',ISYM,
1250     *   ' ROOT NUMBER',IROOT
1251      END IF
1252      RETURN
1253      END
1254C  /* Deck hypind */
1255      SUBROUTINE HYPIND
1256C
1257C CALCULATE A POINTER TO THE NUMBER OF DIFFERENT
1258C LINEAR RESPONSE EQUATIONS THAT NEED TO BE SOLVED
1259C IN A HYPERPOLARIZABILITY CALCULATION
1260C
1261#include "implicit.h"
1262C
1263#include "rspprp.h"
1264#include "infhyp.h"
1265#include "inforb.h"
1266#include "infrsp.h"
1267#include "infpri.h"
1268#include "infspi.h"
1269C
1270      DO 300 ICSYM = 1,NSYM
1271         DO 400 IBSYM = 1,NSYM
1272            IASYM = MULD2H(ICSYM,IBSYM)
1273            IF ( (NCQROP(ICSYM) .GT. 0) .AND. (NBQROP(IBSYM) .GT. 0)
1274     *          .AND. (NAQROP(IASYM).GT.0) ) THEN
1275               DO 500 ICOP = 1,NCQROP(ICSYM)
1276                  DO 600 ICFR = 1,NCQRFR
1277                     IF ( ISPINC.EQ.1 ) THEN
1278                        INUM = INQRTR(CQRLB(ICSYM,ICOP),
1279     *                                CQRFR(ICFR),ICSYM)
1280                     ELSE
1281                        INUM = INQRLR(CQRLB(ICSYM,ICOP),
1282     *                                CQRFR(ICFR),ICSYM)
1283                     END IF
1284 600              CONTINUE
1285 500           CONTINUE
1286               DO 700 IBOP = 1,NBQROP(IBSYM)
1287                  DO 800 IBFR = 1,NBQRFR
1288                     IF ( ISPINB.EQ.1 ) THEN
1289                        INUM = INQRTR(BQRLB(IBSYM,IBOP),
1290     *                                BQRFR(IBFR),IBSYM)
1291                     ELSE
1292                        INUM = INQRLR(BQRLB(IBSYM,IBOP),
1293     *                                BQRFR(IBFR),IBSYM)
1294                     END IF
1295 800              CONTINUE
1296 700           CONTINUE
1297C
1298C  If a special process is specified we only need some of the
1299C  combinations of the frequencies on the B and C operators
1300C
1301               IF (QRSPEC) THEN
1302                  DO IBFR=1,NBQRFR
1303C
1304C  If Pockels effect we need to solve linear response equations for w.
1305C
1306                     IF (QRPOCK) THEN
1307                        DO IAOP=1,NAQROP(IASYM)
1308                           INUM = INQRLR(AQRLB(IASYM,IAOP),
1309     &                                   BQRFR(IBFR),IASYM)
1310                        ENDDO
1311                     ENDIF
1312C
1313C  If SHG we need to solve linear response equations for 2w.
1314C
1315                     IF (QRSHG) THEN
1316                        DO IAOP=1,NAQROP(IASYM)
1317                           INUM = INQRLR(AQRLB(IASYM,IAOP),
1318     &                                 2*BQRFR(IBFR),IASYM)
1319                        ENDDO
1320                       END IF
1321                  ENDDO
1322C
1323C  If optical refractivity we need to solve linear response equations for 0.
1324C
1325                     IF (QROPRF) THEN
1326                        DO IAOP=1,NAQROP(IASYM)
1327                           INUM = INQRLR(AQRLB(IASYM,IAOP),
1328     &                                   0D0,IASYM)
1329                        ENDDO
1330                     ENDIF
1331               ELSE
1332                  DO 1000 ICFR = 1,NCQRFR
1333                     DO 1100 IBFR = 1,NBQRFR
1334                        BPCFR = CQRFR(ICFR) + BQRFR(IBFR)
1335                        DO 1200 IAOP=1,NAQROP(IASYM)
1336                           IF ( ISPINA.EQ.1 ) THEN
1337                              INUM = INQRTR(AQRLB(IASYM,IAOP),
1338     &                                      BPCFR,IASYM)
1339                           ELSE
1340                              INUM = INQRLR(AQRLB(IASYM,IAOP),
1341     &                                      BPCFR,IASYM)
1342                           END IF
1343 1200                   CONTINUE
1344 1100                CONTINUE
1345 1000             CONTINUE
1346               END IF
1347            END IF
1348 400     CONTINUE
1349 300  CONTINUE
1350      RETURN
1351      END
1352C  /* Deck somind */
1353      SUBROUTINE SOMIND
1354C
1355C CALCULATE A POINTER TO THE NUMBER OF DIFFERENT
1356C LINEAR RESPONSE EQUATIONS THAT ARE USED IN A
1357C A CALCULATION OF SECOND ORDER TRANSITION MOMENTS
1358C
1359#include "implicit.h"
1360C
1361#include "rspprp.h"
1362#include "infsmo.h"
1363#include "indqr.h"
1364#include "inforb.h"
1365#include "infrsp.h"
1366#include "infpri.h"
1367#include "infspi.h"
1368C
1369C If two-photon calculation and .BFREQ is not specified
1370C we only compute a restricted number of response vectors
1371C corresponding to BFREQ=CFREQ=EXCITA/2.
1372C
1373      IF (TWOPHO) THEN
1374         NBSMFR = 0
1375         DO ICSYM = 1,NSYM
1376         ICEXC = NSMCNV(ICSYM)
1377         DO IBSYM = 1,NSYM
1378            IASYM = MULD2H(ICSYM,IBSYM)
1379            IF ( (NSMCNV(ICSYM) .GT. 0) .AND. (NBSMOP(IBSYM) .GT. 0)
1380     *          .AND. (NASMOP(IASYM).GT.0) ) THEN
1381               DO IC = 1,NSMCNV(ICSYM)
1382                  IBFR = NBSMFR+IC
1383                  BSMFR(IBFR) = EXCITA(ICSYM,IC,ISPINC+1)/2
1384                  DO IAOP=1,NASMOP(IASYM)
1385                     INUM=INQRLR(ASMLB(IASYM,IAOP),BSMFR(IBFR),IASYM)
1386                  END DO
1387                  DO IBOP=1,NBSMOP(IBSYM)
1388C
1389C It is the response vector Nb(-wf/2) which is needed, but since if
1390C
1391C     Nb(wf/2)  = ( Z Y* ) then Nb(-wf/2) = ( Y* Z )
1392C
1393C so we only compute linear response vectors with positive frequencies.
1394C
1395                     INUM=INQRLR(BSMLB(IBSYM,IBOP),BSMFR(IBFR),IBSYM)
1396                  END DO
1397               END DO
1398               NBSMFR = NBSMFR + ICEXC
1399               ICEXC = 0
1400            END IF
1401         END DO
1402         END DO
1403      ELSE
1404C
1405C The general case computes all combinations.
1406C
1407      DO 300 ICSYM = 1,NSYM
1408         DO 400 IBSYM = 1,NSYM
1409            IASYM = MULD2H(ICSYM,IBSYM)
1410            IF ( (NSMCNV(ICSYM) .GT. 0) .AND. (NBSMOP(IBSYM) .GT. 0)
1411     *          .AND. (NASMOP(IASYM).GT.0) ) THEN
1412               DO 700 IBOP = 1,NBSMOP(IBSYM)
1413                  DO 800 IBFR = 1,NBSMFR
1414                     IF ( ISPINB.EQ.1 ) THEN
1415                        INUM = INQRTR(BSMLB(IBSYM,IBOP),
1416     *                                -BSMFR(IBFR),IBSYM)
1417                     ELSE
1418                        INUM = INQRLR(BSMLB(IBSYM,IBOP),
1419     *                                -BSMFR(IBFR),IBSYM)
1420                     END IF
1421 800              CONTINUE
1422 700           CONTINUE
1423               DO 900 IC = 1,NSMCNV(ICSYM)
1424                  DO 1000 IBFR = 1,NBSMFR
1425                     CEXMBF = EXCITA(ICSYM,IC,ISPINC+1) -BSMFR(IBFR)
1426                     DO 1100 IAOP = 1,NASMOP(IASYM)
1427                        IF ( ISPINA.EQ.0 ) THEN
1428                           INUM = INQRLR(ASMLB(IASYM,IAOP),CEXMBF,IASYM)
1429                        ELSE
1430                           INUM = INQRTR(ASMLB(IASYM,IAOP),CEXMBF,IASYM)
1431                        END IF
1432 1100                CONTINUE
1433 1000             CONTINUE
1434 900           CONTINUE
1435            END IF
1436 400     CONTINUE
1437 300  CONTINUE
1438C
1439      END IF
1440C
1441      RETURN
1442      END
1443C  /* Deck exmind */
1444      SUBROUTINE EXMIND
1445C
1446C CALCULATE A POINTER TO THE NUMBER OF DIFFERENT LINEAR
1447C RESPONSE EQUATIONS THAT ARE USED IN A CALCULATION OF
1448C TRANSITION MOMENTS BETWEEN EXCITED STATES
1449C
1450#include "implicit.h"
1451C
1452#include "infrsp.h"
1453#include "rspprp.h"
1454#include "infpp.h"
1455#include "infsmo.h"
1456#include "indqr.h"
1457#include "inforb.h"
1458#include "infpri.h"
1459#include "infspi.h"
1460#include "priunit.h"
1461C
1462C panor:
1463C     If DOESA=.TRUE, then:
1464C     For excited state absorption we require the B-state to be the
1465C     first state of symmetry ESASYM and that the C-state is separate
1466C     from the B-state.
1467C
1468      DO 300 ICSYM = 1,NSYM
1469         IF (DOESA .OR. EXMTES .OR. ISPINB.NE.ISPINC) THEN
1470            NBSYME = NSYM
1471         ELSE
1472            NBSYME = ICSYM
1473         END IF
1474         DO 400 IBSYM = 1,NBSYME
1475            IF (DOESA.AND.IBSYM.NE.ESASYM) GOTO 400
1476            IASYM = MULD2H(ICSYM,IBSYM)
1477            IF ( (NPPCNV(ICSYM) .GT. 0) .AND. (NPPCNV(IBSYM) .GT. 0)
1478     *          .AND. (NGPPP(IASYM).GT.0) ) THEN
1479               DO 900 IC = 1,NPPCNV(ICSYM)
1480                  IF (DOESA .AND. ICSYM.EQ.IBSYM .AND. IC.EQ.1) GOTO 900
1481                  IF (DOESA) THEN
1482                     IBE = 1
1483                  ELSE IF ( EXMTES .OR. ICSYM.NE.IBSYM .OR.
1484     *               ISPINB .NE. ISPINC) THEN
1485                     IBE = NPPCNV(IBSYM)
1486                  ELSE
1487                     IBE = IC
1488                  END IF
1489                  DO 1000 IB = 1,IBE
1490                     CEXMBE = EXCITA(ICSYM,IC,ISPINC+1) -
1491     *                        EXCITA(IBSYM,IB,ISPINB+1)
1492                     DO 1100 IAOP = 1,NGPPP(IASYM)
1493                        IF (.NOT.RSPCI.AND.
1494     *                       LCMEXC(IB,IBSYM).AND.LCMEXC(IC,ICSYM)) THEN
1495                           IF (ISPINA .EQ. 0) THEN
1496                              INUM = INQRLR
1497     *                               (LBLPP(IASYM,IAOP),CEXMBE,IASYM)
1498                           ELSE
1499                              INUM = INQRTR
1500     *                               (LBLPP(IASYM,IAOP),CEXMBE,IASYM)
1501                           END IF
1502                        END IF
1503 1100                CONTINUE
1504 1000             CONTINUE
1505 900           CONTINUE
1506            END IF
1507 400     CONTINUE
1508 300  CONTINUE
1509      RETURN
1510      END
1511C  /* Deck inqrlr */
1512      INTEGER FUNCTION INQRLR(NEWLBL,FRQVAL,ISYM)
1513C
1514C DETERMINE IF A POINTER TO LINEAR RESPONSE EQUATION
1515C CHARACTERIZED BY NEWLBL AND FRQVAL HAS OCCURRED
1516C IF SO INDQR GIVES ITS NUMBER IN THE LIST
1517C IF THE POINTER DOES NOT APPEAR IT IS ADDED TO
1518C LIST AND INDQR GIVES THE ADDED NUMBER
1519C
1520#include "implicit.h"
1521C
1522#include "priunit.h"
1523#include "infrsp.h"
1524#include "rspprp.h"
1525#include "indqr.h"
1526C
1527      PARAMETER ( DIFFR = 1.0D-8 )
1528      CHARACTER*8 NEWLBL
1529      DO 100 I=1,NLRLBL
1530         IF ( (NEWLBL .EQ. QRLBL(I) ) .AND.
1531     *       ( ABS(FRQVAL- QRFREQ(I)).LE.DIFFR) ) THEN
1532            INQRLR = I
1533            RETURN
1534         END IF
1535 100  CONTINUE
1536      NLRLBL = NLRLBL + 1
1537      IF (NLRLBL.GT.MXLRQR) THEN
1538         WRITE(LUPRI,'(A,/A,I5,A,I5)')
1539     *   ' NUMBER OF SPECIFIED PROPERTIES EXCEED THE MAXIMUM ALLOWED',
1540     *   ' MXLRQR =',MXLRQR,' NLRLBL = ',NLRLBL
1541         CALL QUIT(' INQRLR: TOO MANY PROPERTIES SPECIFIED')
1542      END IF
1543      QRLBL(NLRLBL) = NEWLBL
1544      QRFREQ(NLRLBL) = FRQVAL
1545      ISYMQR(NLRLBL) = ISYM
1546      INQRLR = NLRLBL
1547      IF (IPRRSP.GT.40) THEN
1548         WRITE(LUPRI,'(/A)')' INQRLR: NEW OPERATOR'
1549         WRITE(LUPRI,'(A,I5,2X,A,A,2X,/A,D13.6,2X,A,I5)')
1550     *   ' NUMBER:',NLRLBL,' LABEL: ',NEWLBL,' FREQUENCY: ',FRQVAL,
1551     *   'SYMMETRY',ISYM
1552      END IF
1553      RETURN
1554      END
1555C  /* Deck qrlrsr */
1556      SUBROUTINE QRLRSR
1557C
1558C SORT THE LINEAR RESPONSE EQUATIONS AFTER SYMMETRY
1559C THE NUMBER OF LINEAR EQUATIONS OF A GIVEN SYMMETRY ARE
1560C STORED IN NLRQR WITH OFSET IN ILRQR
1561C
1562#include "implicit.h"
1563C
1564      CHARACTER*8 LABEL
1565C
1566#include "priunit.h"
1567#include "inforb.h"
1568#include "rspprp.h"
1569#include "indqr.h"
1570#include "infpri.h"
1571#include "infrsp.h"
1572C
1573      IF (IPRRSP.GT.40) THEN
1574         WRITE(LUPRI,'(/A,I5,2X,A)')
1575     *   ' LIST OF NLRLBL= ',NLRLBL,
1576     *   ' LINEAR EQUATIONS FOR QUADRATIC RESPONSE'
1577         WRITE(LUPRI,'(A)')' BEFORE SORTING'
1578         WRITE(LUPRI,'(/A)')'   I  QRLBL(I)  QRFREQ(I)  ISYMQR(I) '
1579         DO 90 I = 1,NLRLBL
1580            WRITE(LUPRI,'(I5,5X,A,D13.6,I12)')
1581     *      I,QRLBL(I),QRFREQ(I),ISYMQR(I)
1582  90     CONTINUE
1583      END IF
1584      ISMMAX = 1
1585      DO ISYM = 1,NSYM
1586         DO ISRT = 1, NLRLBL
1587            IF (ISRT.GE.ISMMAX .AND. ISYMQR(ISRT).EQ.ISYM) THEN
1588               IQR   = ISYMQR(ISRT)
1589               FRVAL = QRFREQ(ISRT)
1590               LABEL = QRLBL(ISRT)
1591               ISYMQR(ISRT) = ISYMQR(ISMMAX)
1592               QRFREQ(ISRT) = QRFREQ(ISMMAX)
1593               QRLBL(ISRT)  = QRLBL(ISMMAX)
1594               ISYMQR(ISMMAX) = IQR
1595               QRFREQ(ISMMAX) = FRVAL
1596               QRLBL(ISMMAX)  = LABEL
1597               ISMMAX = ISMMAX + 1
1598            END IF
1599         END DO
1600      END DO
1601      ITOTOP = 0
1602      DO 300 ISYM = 1,NSYM
1603         INUMOP = 0
1604         DO 400 I = 1,NLRLBL
1605            IF(ISYMQR(I).EQ.ISYM) THEN
1606               INUMOP = INUMOP + 1
1607            END IF
1608 400     CONTINUE
1609         ILRQR(ISYM) = ITOTOP
1610         ITOTOP = ITOTOP + INUMOP
1611         NLRQR(ISYM) = INUMOP
1612 300  CONTINUE
1613      IF (IPRRSP.GT.40) THEN
1614         WRITE(LUPRI,'(/A,I5,2X,A)')
1615     *   ' LIST OF NLRLBL= ',NLRLBL,
1616     *   ' LINEAR EQUATIONS FOR QUADRATIC RESPONSE'
1617         WRITE(LUPRI,'(A)')'AFTER SORTING'
1618         DO 301 I = 1,NLRLBL
1619            WRITE(LUPRI,'(/A)')'   I  QRLBL(I)  QRFREQ(I)  ISYMQR(I) '
1620            WRITE(LUPRI,'(I5,5X,A,D13.6,I12)')
1621     *      I,QRLBL(I),QRFREQ(I),ISYMQR(I)
1622 301     CONTINUE
1623      END IF
1624      IF (IPRRSP.GT.35) THEN
1625         WRITE(LUPRI,'(/A)')
1626     *   ' NUMBER OF LINEAR RESPONSE EQUATIONS IN VARIOUS SYMMETRIES'
1627         WRITE(LUPRI,'(/A)')' NLRQR(ISYM),ISYM=1,NSYM '
1628         WRITE(LUPRI,'(8I5)')( NLRQR(ISYM),ISYM=1,NSYM )
1629      END IF
1630      RETURN
1631      END
1632C  /* Deck qrexsr */
1633      SUBROUTINE QREXSR
1634C
1635C SORT THE EXCITATION ENERGIES  AFTER SYMMETRY
1636C THE NUMBER OF EXCITATION ENERGIES OF A GIVEN SYMMETRY ARE
1637C STORED IN NEXQR WITH OFFSET IN IEXQR
1638C
1639#include "implicit.h"
1640C
1641#include "priunit.h"
1642#include "infrsp.h"
1643#include "inforb.h"
1644#include "rspprp.h"
1645#include "indqr.h"
1646#include "infpri.h"
1647C
1648      IF (IPRRSP.GT.40) THEN
1649         WRITE(LUPRI,'(/A,I5,2X,A)')
1650     *   ' LIST OF NEXLBL= ',NEXLBL,'EXCITATIONS FOR QUADRATIC RESPONSE'
1651         WRITE(LUPRI,'(A)')' BEFORE SORTING'
1652         DO 90 I = 1,NEXLBL
1653            WRITE(LUPRI,'(/A)')'     I   ISEXQR(I)  JEXQR(I)'
1654            WRITE(LUPRI,'(/3I10)')I,ISEXQR(I),JEXQR(I)
1655 90      CONTINUE
1656      END IF
1657      ISMMAX = 1
1658      DO 100 ISYM = 1,NSYM
1659         I = ISMMAX
1660         DO 200 ISRT = I,NEXLBL
1661            IF ( ISEXQR(ISRT).EQ.ISYM) THEN
1662                IQR   = ISEXQR(ISRT)
1663                JQR   = JEXQR(ISRT)
1664                ISEXQR(ISRT) = ISEXQR(ISMMAX)
1665                JEXQR(ISRT)  = JEXQR(ISMMAX)
1666                ISEXQR(ISMMAX) = IQR
1667                JEXQR(ISMMAX)  = JQR
1668                ISMMAX = ISMMAX + 1
1669            END IF
1670 200     CONTINUE
1671 100  CONTINUE
1672      ITOTOP = 0
1673      DO 300 ISYM = 1,NSYM
1674         INUMOP = 0
1675         DO 400 I = 1,NEXLBL
1676            IF(ISEXQR(I).EQ.ISYM) THEN
1677               INUMOP = INUMOP + 1
1678            END IF
1679 400     CONTINUE
1680         IEXQR(ISYM) = ITOTOP
1681         ITOTOP = ITOTOP + INUMOP
1682         NEXQR(ISYM) = INUMOP
1683 300  CONTINUE
1684      IF (IPRRSP.GT.40) THEN
1685         WRITE(LUPRI,'(/A,I5,2X,A)')
1686     *   ' LIST OF NEXLBL= ',NEXLBL,'EXCITATIONS FOR QUADRATIC RESPONSE'
1687         WRITE(LUPRI,'(A)')' AFTER SORTING'
1688         DO 310 I = 1,NEXLBL
1689            WRITE(LUPRI,'(/A)')'     I   ISEXQR(I)  JEXQR(I)'
1690            WRITE(LUPRI,'(/3I10)')I,ISEXQR(I),JEXQR(I)
1691 310     CONTINUE
1692      END IF
1693      IF (IPRRSP.GT.35) THEN
1694         WRITE(LUPRI,'(/A)')
1695     *   ' NUMBER OF EXCITATIONS IN VARIOUS SYMMETRIES'
1696         WRITE(LUPRI,'(/A)')' NEXQR(ISYM),ISYM=1,NSYM '
1697         WRITE(LUPRI,'(8I5)')( NEXQR(ISYM),ISYM=1,NSYM )
1698      END IF
1699      RETURN
1700      END
1701C  /* Deck inqrtr */
1702      INTEGER FUNCTION INQRTR(NEWLBL,FRQVAL,ISYM)
1703C
1704C DETERMINE IF A POINTER TO TRIPLET LINEAR RESPONSE EQUATION
1705C CHARACTERIZED BY NEWLBL AND FRQVAL
1706C HAS OCCURED IF SO INQRTR GIVES ITS NUMBER IN THE
1707C LIST IF THE POINTER DOES NOT APPEAR IT IS ADDED TO
1708C LIST AND INQRTR GIVES THE ADDED NUMBER
1709C
1710#include "implicit.h"
1711C
1712#include "priunit.h"
1713#include "infrsp.h"
1714#include "rspprp.h"
1715#include "indqr.h"
1716#include "infpri.h"
1717C
1718      PARAMETER ( DIFFR = 1.0D-8 )
1719      CHARACTER*8 NEWLBL
1720      DO 100 I=1,NTRLBL
1721         IF ( (NEWLBL .EQ. TRLBL(I) ) .AND.
1722     *       ( ABS(FRQVAL- TRFREQ(I)).LE.DIFFR) ) THEN
1723            INQRTR = I
1724            RETURN
1725         END IF
1726 100  CONTINUE
1727      NTRLBL = NTRLBL + 1
1728      IF (NTRLBL.GT.MXLRQR) THEN
1729         WRITE(LUPRI,'(A,/A,I5,A,I5)')
1730     *   ' NUMBER OF SPECIFIED PROPERTIES EXCEED THE MAXIMUM ALLOWED',
1731     *   ' MXLRQR =',MXLRQR,' NTRLBL = ',NTRLBL
1732         CALL QUIT(' INQRTR: TOO MANY PROPERTIES SPECIFIED')
1733      END IF
1734      TRLBL(NTRLBL) = NEWLBL
1735      TRFREQ(NTRLBL) = FRQVAL
1736      ISYMTR(NTRLBL) = ISYM
1737      INQRTR = NTRLBL
1738      IF (IPRRSP.GT.40) THEN
1739         WRITE(LUPRI,'(/A)')' INQRTR: NEW OPERATOR'
1740         WRITE(LUPRI,'(A,I5,2X,A,A,2X,/A,D13.6,2X,A,I5)')
1741     *   ' NUMBER:',NTRLBL,' LABEL: ',NEWLBL,' FREQUENCY: ',FRQVAL,
1742     *   'SYMMETRY',ISYM
1743      END IF
1744      RETURN
1745      END
1746C  /* Deck qrtrsr */
1747      SUBROUTINE QRTRSR
1748C
1749C SORT THE TRIPLET LINEAR RESPONSE EQUATIONS AFTER SYMMETRY
1750C THE NUMBER OF LINEAR EQUATIONS OF A GIVEN SYMMETRY ARE
1751C STORED IN NTRQR WITH OFSET IN ITRQR
1752C
1753#include "implicit.h"
1754C
1755      CHARACTER*8 LABEL
1756C
1757#include "priunit.h"
1758#include "inforb.h"
1759#include "rspprp.h"
1760#include "indqr.h"
1761#include "infpri.h"
1762#include "infrsp.h"
1763C
1764      IF (IPRRSP.GT.40) THEN
1765         WRITE(LUPRI,'(/A,I5,2X,A)')
1766     *   ' LIST OF NTRLBL= ',NTRLBL,
1767     *   ' TRIPLET LINEAR EQUATIONS FOR QUADRATIC RESPONSE'
1768         WRITE(LUPRI,'(A)')' BEFORE SORTING'
1769         WRITE(LUPRI,'(/A)')'   I  TRLBL(I)  TRFREQ(I)  ISYMTR(I) '
1770         DO 90 I = 1,NTRLBL
1771            WRITE(LUPRI,'(I5,5X,A,D13.6,I12)')
1772     *      I,TRLBL(I),TRFREQ(I),ISYMTR(I)
1773  90     CONTINUE
1774      END IF
1775      ISMMAX = 1
1776      DO 100 ISYM = 1,NSYM
1777         I = ISMMAX
1778         DO 200 ISRT = I,NTRLBL
1779            IF ( ISYMTR(ISRT).EQ.ISYM) THEN
1780                IQR   = ISYMTR(ISRT)
1781                FRVAL  = TRFREQ(ISRT)
1782                LABEL = TRLBL(ISRT)
1783                ISYMTR(ISRT) = ISYMTR(ISMMAX)
1784                TRFREQ(ISRT) = TRFREQ(ISMMAX)
1785                TRLBL(ISRT)  = TRLBL(ISMMAX)
1786                ISYMTR(ISMMAX) = IQR
1787                TRFREQ(ISMMAX) = FRVAL
1788                TRLBL(ISMMAX)  = LABEL
1789                ISMMAX = ISMMAX + 1
1790            END IF
1791 200     CONTINUE
1792 100  CONTINUE
1793      ITOTOP = 0
1794      DO 300 ISYM = 1,NSYM
1795         INUMOP = 0
1796         DO 400 I = 1,NTRLBL
1797            IF(ISYMTR(I).EQ.ISYM) THEN
1798               INUMOP = INUMOP + 1
1799            END IF
1800 400     CONTINUE
1801         ITRQR(ISYM) = ITOTOP
1802         ITOTOP = ITOTOP + INUMOP
1803         NTRQR(ISYM) = INUMOP
1804 300  CONTINUE
1805      IF (IPRRSP.GT.40) THEN
1806         WRITE(LUPRI,'(/A,I5,2X,A/A)')
1807     *   ' LIST OF NTRLBL= ',NTRLBL,
1808     *   ' TRIPLET LINEAR EQUATIONS FOR QUADRATIC RESPONSE',
1809     *   ' AFTER SORTING'
1810         DO 301 I = 1,NTRLBL
1811            WRITE(LUPRI,'(/A)')'   I  TRLBL(I)  TRFREQ(I)  ISYMTR(I) '
1812            WRITE(LUPRI,'(I5,5X,A,D13.6,I12)')
1813     *      I,TRLBL(I),TRFREQ(I),ISYMTR(I)
1814 301     CONTINUE
1815      END IF
1816      IF (IPRRSP.GT.35) THEN
1817         WRITE(LUPRI,'(/2A)')
1818     *   ' NUMBER OF TRIPLET LINEAR RESPONSE EQUATIONS',
1819     *   ' IN VARIOUS SYMMETRIES'
1820         WRITE(LUPRI,'(/A)')' NTRQR(ISYM),ISYM=1,NSYM '
1821         WRITE(LUPRI,'(8I5)')( NTRQR(ISYM),ISYM=1,NSYM )
1822      END IF
1823      RETURN
1824      END
1825