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
19*=====================================================================*
20      SUBROUTINE CC_SOLRSP(WORK,LWORK,APROXR12)
21C---------------------------------------------------------------------*
22C
23C     Purpose: call equation driver for the various types of
24C              response equations we have to solve
25C
26C      Written by Ove Christiansen October 1996
27C      Restructured by Christof Haettig, 1997/1998
28C      PL1. Sonia 2000
29C      SLV98,OC: Self consistent solution of 1-order response equations.
30C
31C---------------------------------------------------------------------*
32C
33      USE PELIB_INTERFACE, ONLY: USE_PELIB
34      IMPLICIT NONE
35#include "priunit.h"
36#include "maxorb.h"
37#include "dummy.h"
38#include "iratdef.h"
39#include "cclr.h"
40#include "ccorb.h"
41#include "ccsdsym.h"
42#include "ccsdio.h"
43#include "ccsdinp.h"
44#include "ccsections.h"
45#include "ccexci.h"
46#include "cclrmrsp.h"
47#include "ccer1rsp.h"
48#include "ccer2rsp.h"
49#include "ccel1rsp.h"
50#include "ccel2rsp.h"
51#include "ccr1rsp.h"
52#include "cco1rsp.h"
53#include "ccx1rsp.h"
54#include "ccl1rsp.h"
55#include "ccr2rsp.h"
56#include "ccl2rsp.h"
57#include "cco2rsp.h"
58#include "ccx2rsp.h"
59#include "ccr3rsp.h"
60#include "ccl3rsp.h"
61#include "cco3rsp.h"
62#include "ccx3rsp.h"
63#include "ccr4rsp.h"
64#include "ccl4rsp.h"
65#include "cco4rsp.h"
66#include "ccx4rsp.h"
67#include "ccn2rsp.h"
68#include "ccrc1rsp.h"
69#include "cclc1rsp.h"
70#include "cccr2rsp.h"
71#include "cccl2rsp.h"
72#include "ccco2rsp.h"
73#include "cccx2rsp.h"
74#include "ccexgr.h"
75#include "leinf.h"
76#include "ccpl1rsp.h"
77#include "ccslvinf.h"
78#include "cctpainf.h"
79Cholesky
80#include "ccdeco.h"
81Cholesky
82C
83C
84      LOGICAL LOCDBG, lold
85      PARAMETER (LOCDBG = .FALSE., lold = .false.)
86C
87      CHARACTER*(*) APROXR12
88      CHARACTER*8 FC1AM,FC2AM,FRHO1,FRHO2
89      PARAMETER (FC1AM='CCR_C1AM',FC2AM='CCR_C2AM')
90      PARAMETER (FRHO1='CCR_RHO1',FRHO2='CCR_RHO2')
91
92      INTEGER LWORK
93
94#if defined (SYS_CRAY)
95      REAL TOL, ONE, TWO, WORK(LWORK), DDOT, X1, X2, FREK
96#else
97      DOUBLE PRECISION TOL, ONE, TWO, WORK(LWORK), DDOT, X1, X2, FREK
98#endif
99      LOGICAL LDUMMY
100      INTEGER ISYM0, IZERO
101
102      PARAMETER(TOL=1.0D-12, ISYM0 = 1, LDUMMY = .FALSE.)
103      PARAMETER(ONE = 1.0D0, IZERO = 0, TWO = 2.0D0)
104
105      CHARACTER MODEL*10,MODEL1*10
106      CHARACTER*3 LIST
107      INTEGER MAXFDIM
108      PARAMETER (MAXFDIM = 500)
109      INTEGER IFTRAN(3,MAXFDIM)
110      INTEGER ORDER, I1DXORD, IOPREL, NSTAT, ISIDE, IDXR1, IOPT,
111     &        ILSTNR, ISYM, KT1, KT2, KEND1, LEND1, IDXM, IDXRE, IDXF,
112     &        NTAMP, KETA, LWRK1, KBEVC, KEND2, LWRK2, IDXRC, KVEC1,
113     &        KVEC2, NCAU, IVEC, INUM, IDXLC, IDXL2, IDXR2, IDXL3,
114     &        IDXR3, NLRCLBLSAVE, KTR12
115C
116      LOGICAL R1HFSKIP, O1SKIP, X1SKIP
117      SAVE    R1HFSKIP, O1SKIP, X1SKIP
118      DATA R1HFSKIP / .FALSE. /
119      DATA O1SKIP   / .FALSE. /
120      DATA X1SKIP   / .FALSE. /
121
122* external functions:
123      INTEGER IR1TAMP, ILRCAMP, ILC1AMP, ICR2AMP, IRHSCR2, IR2TAMP,
124     &        ICL2AMP, IETACL2, IR3TAMP, ILSTSYM
125
126C-------------------------------------
127C     Some initializations:
128C-------------------------------------
129C
130      CALL QENTER('CC_SOLRSP')
131C
132      MODEL = 'CCSD      '
133      IF (CCS) MODEL = 'CCS       '
134      IF (CC2) MODEL = 'CC2       '
135C
136      IF (.NOT.(CCS.OR.CC2.OR.CCSD.OR.CC3)) THEN
137         WRITE(LUPRI,'(A)') ' CC_SOLRSP: Do not want to calculate '
138     *           //'anything else than CCS, CC2 and CCSD properties '
139         CALL QUIT('Model not CCS or CC2 or CCSD in CC_SOLRSP')
140      ENDIF
141
142      IF (CHOINT .AND. (.NOT. CC2)) THEN
143         WRITE(LUPRI,'(A,A,A)')
144     &   ' CC_SOLRSP: Cholesky decomposed integrals not implemented',
145     &   ' for ',MODEL
146         CALL QUIT('Cholesky '//MODEL//' not implemented in CC_SOLRSP')
147      ENDIF
148
149      WRITE (LUPRI,'(4X,A,/)') '  '
150      WRITE (LUPRI,'(4X,A)')
151     *'*********************************************************'//
152     *'**********'
153      WRITE (LUPRI,'(4X,A)')
154     *'*          SOLVING COUPLED CLUSTER RESPONSE EQUATIONS    '//
155     *'         *'
156      WRITE (LUPRI,'(4X,A,///)')
157     *'*********************************************************'//
158     *'**********'
159
160      IF (CHOINT) THEN
161         WRITE(LUPRI,'(4X,A,///)')
162     &   '         (Calculation using Cholesky decomposed integrals)'
163      ELSE
164         WRITE(LUPRI,'(A)') ' '
165      ENDIF
166
167C======================================================================
168
169      WRITE(LUPRI,'(1X,A1,70("="),A1)') '+','+'
170      WRITE(LUPRI,'(1X,A1,21X,A,20X,A1)')
171     &  '!','RHS & ETA VECTORS TO COMPUTE:','!'
172      WRITE(LUPRI,'(1X,A1,70("="),A1)') '+','+'
173
174      WRITE(LUPRI, '(1X,A1,1X,A4,1X,A1,1X,A6,1X,A1,1X,A52,1X,A1)')
175     &   '|','TYPE','|','# VEC.','|',
176     &   'NEEDED FOR:                                        ','|'
177
178      WRITE(LUPRI,'(1X,A1,70("-"),A1)') '+','+'
179
180      IF (NO1LBL.GT.0) WRITE(LUPRI,2) '|',' O1 ','|',   NO1LBL,'|',
181     &   'first-order amplitude equations                    ','|'
182      IF (NX1LBL.GT.0) WRITE(LUPRI,2) '|',' X1 ','|',   NX1LBL,'|',
183     &   'first-order multipliers equations                  ','|'
184      IF (NO2LBL.GT.0) WRITE(LUPRI,2) '|',' O2 ','|',   NO2LBL,'|',
185     &   'second-order amplitude equations                   ','|'
186      IF (NX2LBL.GT.0) WRITE(LUPRI,2) '|',' X2 ','|',   NX2LBL,'|',
187     &   'second-order multipliers equations                 ','|'
188      IF (NCO2LBL.GT.0) WRITE(LUPRI,2) '|','CO2','|',  NCO2LBL,'|',
189     &   'second-order right Cauchy equations                ','|'
190      IF (NCX2LBL.GT.0) WRITE(LUPRI,2) '|','CX2','|',  NCX2LBL,'|',
191     &   'second-order left Cauchy equations                 ','|'
192      IF (NO3LBL.GT.0) WRITE(LUPRI,2) '|',' O3 ','|',   NO3LBL,'|',
193     &   'third-order amplitude equations                    ','|'
194      IF (NX3LBL.GT.0) WRITE(LUPRI,2) '|',' X3 ','|',   NX3LBL,'|',
195     &   'third-order multiplier equations                   ','|'
196      IF (NXGRST.GT.0) WRITE(LUPRI,2) '|',' BE ','|',   NXGRST,'|',
197     &   'zeroth-order excited-state Lagrange multiplier     ','|'
198      IF (NQRN2.GT.0)  WRITE(LUPRI,2) '|',' BR ','|',    NQRN2,'|',
199     &   'zeroth-order excited-state transition Lagr. multip.','|'
200      IF (NER1LBL.GT.0) WRITE(LUPRI,2)'|',' EO1','|',  NER1LBL,'|',
201     &   'first-order excited state right response equations ','|'
202      IF (NEL1LBL.GT.0) WRITE(LUPRI,2)'|',' EX1','|',  NEL1LBL,'|',
203     &   'first-order excited state left response equations  ','|'
204      IF (NEL2LBL.GT.0) WRITE(LUPRI,2)'|',' EO2','|',  NER2LBL,'|',
205     &   'second-order excited state right response equations','|'
206      IF (NER2LBL.GT.0) WRITE(LUPRI,2)'|',' EX2','|',  NEL2LBL,'|',
207     &   'second-order excited state left response equations ','|'
208
209      WRITE(LUPRI,'(1X,A1,70("="),A1,///)') '+','+'
210
211C======================================================================
212
213      WRITE(LUPRI,'(1X,A1,70("="),A1)') '+','+'
214      WRITE(LUPRI,'(1X,A1,21X,A,20X,A1)')
215     &  '!','  LINEAR EQUATIONS TO SOLVE: ','!'
216      WRITE(LUPRI,'(1X,A1,70("="),A1)') '+','+'
217
218      WRITE(LUPRI, '(1X,A1,1X,A4,1X,A1,1X,A6,1X,A1,1X,A52,1X,A1)')
219     &   '|','TYPE','|','# VEC.','|',
220     &   'EQUATION:                                          ','|'
221
222      WRITE(LUPRI,'(1X,A1,70("-"),A1)') '+','+'
223
224      IF (NLRTHFLBL.GT.0)WRITE(LUPRI,2)'|',' R1 ','|',NLRTHFLBL,'|',
225     &   'first-order orbtial response                       ','|'
226      IF (NLRTLBL.GT.0) WRITE(LUPRI,2) '|',' R1 ','|', NLRTLBL,'|',
227     &   'first-order amplitude response                     ','|'
228      IF (NLRZLBL.GT.0) WRITE(LUPRI,2) '|',' L1 ','|', NLRZLBL,'|',
229     &   'first-order multiplier response                    ','|'
230      IF (NLRCLBL.GT.0) WRITE(LUPRI,2) '|',' RC ','|', NLRCLBL,'|',
231     &   'first-order right Cauchy vectors                   ','|'
232      IF (NLC1LBL.GT.0) WRITE(LUPRI,2) '|',' LC ','|', NLC1LBL,'|',
233     &   'first-order left Cauchy vectors                    ','|'
234      IF (NR2TLBL.GT.0) WRITE(LUPRI,2) '|',' R2 ','|', NR2TLBL,'|',
235     &   'second-order amplitude response                    ','|'
236      IF (NCR2LBL.GT.0) WRITE(LUPRI,2) '|',' CR2','|', NCR2LBL,'|',
237     &   'second-order right Cauchy vectors                  ','|'
238      IF (NL2LBL.GT.0)  WRITE(LUPRI,2) '|',' L2 ','|',  NL2LBL,'|',
239     &   'second-order multiplier response                   ','|'
240      IF (NCL2LBL.GT.0) WRITE(LUPRI,2) '|',' CL2','|', NCL2LBL,'|',
241     &   'second-order left Cauchy vectors                   ','|'
242      IF (NR3TLBL.GT.0) WRITE(LUPRI,2) '|',' R3 ','|', NR3TLBL,'|',
243     &   'third-order amplitude response                     ','|'
244      IF (NL3LBL.GT.0)  WRITE(LUPRI,2) '|',' L3 ','|',  NL3LBL,'|',
245     &   'third-order multiplier response                    ','|'
246      IF (NLRM.GT.0)    WRITE(LUPRI,2) '|',' M1 ','|',    NLRM,'|',
247     &   'first-order ground - excited state multipliers     ','|'
248      IF (NXGRST.GT.0)  WRITE(LUPRI,2) '|',' E0 ','|',  NXGRST,'|',
249     &   'zeroth-order excited-state Lagrange multiplier     ','|'
250      IF (NQRN2.GT.0)   WRITE(LUPRI,2) '|',' N2 ','|',   NQRN2,'|',
251     &   'zeroth-order excited-state transition Lagr. multip.','|'
252      IF (NER1LBL.GT.0) WRITE(LUPRI,2) '|',' ER1','|', NER1LBL,'|',
253     &   'first-order right excited state response           ','|'
254      IF (NEL1LBL.GT.0) WRITE(LUPRI,2) '|',' EL1','|', NEL1LBL,'|',
255     &   'first-order left excited state response            ','|'
256      IF (NER2LBL.GT.0) WRITE(LUPRI,2) '|',' ER2','|', NER2LBL,'|',
257     &   'second-order right excited state response          ','|'
258      IF (NEL2LBL.GT.0) WRITE(LUPRI,2) '|',' EL2','|', NEL2LBL,'|',
259     &   'second-order left excited state response           ','|'
260      IF (NPL1LBL.GT.0) WRITE(LUPRI,2) '|',' PL1','|', NPL1LBL,'|',
261     &   'first-order projected multiplier response          ','|'
262
263      WRITE(LUPRI,'(1X,A1,70("="),A1,///)') '+','+'
264
265C======================================================================
266
267      WRITE(LUPRI,'(1X,A1,70("="),A1)') '+','+'
268      WRITE(LUPRI,'(1X,A1,16X,A,16X,A1)')
269     &  '!','F MATRIX TRANSFORMATIONS TO COMPUTE:  ','!'
270      WRITE(LUPRI,'(1X,A1,70("="),A1)') '+','+'
271
272      WRITE(LUPRI, '(1X,A1,1X,A4,1X,A1,1X,A6,1X,A1,1X,A52,1X,A1)')
273     &   '|','TYPE','|','# VEC.','|',
274     &   'TRANSFORMED:                                       ','|'
275
276      WRITE(LUPRI,'(1X,A1,70("-"),A1)') '+','+'
277
278      IF (NLRTLBL.GT.0) WRITE(LUPRI,2) '|',' F1 ','|', NLRTLBL,'|',
279     &   'first-order amplitude response (R1) vector         ','|'
280      IF (NLC1LBL.GT.0) WRITE(LUPRI,2) '|',' FC ','|', NLC1LBL,'|',
281     &   'first-order right Cauchy (RC) vector               ','|'
282      IF (NLRM.GT.0)    WRITE(LUPRI,2) '|',' FR ','|',    NLRM,'|',
283     &   'right eigenvectors (RE)                            ','|'
284      IF (NL2LBL.GT.0)  WRITE(LUPRI,2) '|',' F2 ','|',  NL2LBL,'|',
285     &   'second-order amplitude response (R2) vector        ','|'
286      IF (NCL2LBL.GT.0) WRITE(LUPRI,2) '|',' CF2','|', NCL2LBL,'|',
287     &   'second-order right Cauchy (CR2) vector             ','|'
288      IF (NL3LBL.GT.0)  WRITE(LUPRI,2) '|',' F3 ','|',  NL3LBL,'|',
289     &   'third-order amplitude response (R3) vector         ','|'
290
291      WRITE(LUPRI,'(1X,A1,70("="),A1,///)') '+','+'
292
293C======================================================================
294C     Solve first-order CPHF kappa equations.
295C     (since these equations are independent of the CC amplitudes,
296C      we solve them only the first time we enter this routine)
297C======================================================================
298      IF ( .NOT. R1HFSKIP ) THEN
299
300         CALL CC_CPHF('R1 ',LRTHFLBL,IDUMMY,IDUMMY,DUMMY,
301     &                ISYLRTHF,FRQLRTHF,IDUMMY,NLRTHFLBL,MAXTLBL,
302     &                WORK,LWORK)
303
304         R1HFSKIP = .TRUE.
305
306      END IF
307
308C======================================================================
309C    Calculate effective Fock matrices from Kappa^(1) one-index
310C    transformed integrals:
311C======================================================================
312
313C     IF (NSYM.EQ.1) THEN
314        I1DXORD = 1
315        IOPREL  = 0
316        CALL CC_GRAD2(I1DXORD,WORK,LWORK)
317C       CALL CCEFFFOCK(I1DXORD,IOPREL,WORK,LWORK)
318C     ELSE
319C       WRITE (LUPRI,*) 'Warning: no 1-index transformed eff. Fock '
320C       WRITE (LUPRI,*) 'matrices calculated... '
321C       WRITE (LUPRI,*)
322C    *             'RELAXED second-order properties will be wrong!.'
323C     END IF
324
325C======================================================================
326C     precalculate right-hand-sides vectors "O1" and the eta vectors
327C     for those first-order cluster amplitude and Lagrangian multiplier
328C     equations for which this is neccessary, because
329C     they involve orbital relaxaton and/or derivative integrals
330C     precalculate first-order chi vectors "X1" for those first-order
331C======================================================================
332Ctmp
333      CALL GPINQ('SKIP_O1','EXIST',O1SKIP)
334Ctmp
335      IF ( .NOT. O1SKIP ) THEN
336
337          CALL CCRHSVEC('O1 ',LBLO1,IDUMMY,IDUMMY,DUMMY,ISYO1,FRQO1,
338     &                  LORXO1,IDUMMY,NO1LBL,MAXO1LBL,0,
339     &                  WORK,LWORK)
340
341      ELSE
342        WRITE (LUPRI,*)
343     &        'Skipped calculation of O1 right hand side vectors.'
344        O1SKIP = .FALSE.
345        CALL SYSTEM('rm SKIP_O1')
346      END IF
347C
348C------------------------------------------------------------------
349C     Start loop over self-consistent solvent solution of R1/L1
350C     CCSLV98, OC
351C------------------------------------------------------------------
352C
353      ICCSLIT = 0
354C
355  999 CONTINUE
356C======================================================================
357C     Solve right first-order cluster amplitude response equations.
358C======================================================================
359      LIST  = 'R1 '
360      NSTAT = 0
361      ORDER = 1
362      ISIDE = +1
363
364      IF ( .NOT. R1SKIP ) THEN
365         CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
366     &                  IDUMMY, IDUMMY, DUMMY, LDUMMY,
367     &                  ISYLRT, LRTLBL, FRQLRT, IDUMMY,
368     &                  ISYOFT, NLRTLBL, MAXTLBL,
369     &                  WORK, LWORK)
370
371      ELSE
372         WRITE (LUPRI,*) 'Skipped first-order response (R1) equations'
373         R1SKIP = .FALSE.
374      END IF
375
376      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
377         CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
378     &                  IDUMMY, IDUMMY, DUMMY, LDUMMY,
379     &                  ISYLRT, LRTLBL, FRQLRT, IDUMMY,
380     &                  ISYOFT, NLRTLBL, MAXTLBL, WORK, LWORK)
381      END IF
382
383C======================================================================
384C     precalculate F * R1 transformation required in CCLR section
385C     and to for the first-order coupled cluster multipliers L1
386C======================================================================
387      LIST  = 'F1 '
388      NSTAT = 0
389      ORDER = 1
390      ISIDE = +1
391
392      IF ( .NOT. F1SKIP ) THEN
393
394         IF (NLRTLBL .GT. 0) THEN
395
396            IF (NLRTLBL .GT. MAXFDIM) THEN
397              WRITE (LUPRI,*)
398     &              'ERROR: Dimension MAXFDIM to small in CC_SOLRSP.'
399              CALL QUIT('Dimension MAXFDIM to small in CC_SOLRSP.')
400            END IF
401C
402C           set up list of F matrix transformations:
403C
404            DO IDXR1 = 1, NLRTLBL
405              IFTRAN(1,IDXR1) = 0
406              IFTRAN(2,IDXR1) = IDXR1
407              IFTRAN(3,IDXR1) = IDXR1
408            END DO
409C
410C           call CC_FMATRIX to do the F * R1 transformations:
411C
412            IOPT = 3
413            CALL CC_FMATRIX(IFTRAN, NLRTLBL, 'L0 ', 'R1 ', IOPT, 'F1 ',
414     &                      IDUMMY, DUMMY, 0, WORK, LWORK )
415
416         END IF
417
418      ELSE
419        WRITE (LUPRI,*) 'Skipped of F x R1 transformations'
420        F1SKIP = .FALSE.
421      END IF
422
423
424      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
425         CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
426     &                  IDUMMY, IDUMMY, DUMMY, LDUMMY,
427     &                  ISYLRT, LRTLBL, FRQLRT, IDUMMY,
428     &                  ISYOFT, NLRTLBL, MAXTLBL, WORK, LWORK)
429      END IF
430
431C======================================================================
432C     For Cholesky CC2, calculate the "right" first-order density
433C     matrix.
434C======================================================================
435
436      IF (CHOINT .AND. CC2) THEN
437         IF (.NOT. D01SKIP) THEN
438
439            ORDER = 1
440            ISIDE = 1
441
442            DO IDXR1 = 1, NLRTLBL
443               IFTRAN(1,IDXR1) = 0
444               IFTRAN(2,IDXR1) = IDXR1
445               IFTRAN(3,IDXR1) = IDXR1
446            ENDDO
447
448            CALL CC_CHOPDEN(IFTRAN,NLRTLBL,'L0 ','R1 ','d01',ORDER,
449     &                      ISIDE,'DENSITY_01',.TRUE.,WORK,LWORK)
450
451         ELSE
452
453            WRITE (LUPRI,*) 'Skipped right first-order density (d01)',
454     &                      ' calculation'
455            D01SKIP = .FALSE.
456
457         ENDIF
458      ENDIF
459
460C======================================================================
461C     Solve first-order Lagrange multiplier response equations.
462C======================================================================
463      LIST  = 'L1 '
464      NSTAT = 0
465      ORDER = 1
466      ISIDE = -1
467
468      IF ( .NOT. L1SKIP ) THEN
469        CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
470     &                 IDUMMY, IDUMMY, DUMMY, LDUMMY,
471     &                 ISYLRZ, LRZLBL, FRQLRZ, IDUMMY,
472     &                 ISYOFZ, NLRZLBL, MAXZLBL,
473     &                 WORK, LWORK)
474
475      ELSE
476        WRITE (LUPRI,*)
477     &        'Skipped first-order left response (L1) equations'
478        L1SKIP = .FALSE.
479      END IF
480
481      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
482        CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
483     &                 IDUMMY, IDUMMY, DUMMY, LDUMMY,
484     &                 ISYLRZ, LRZLBL, FRQLRZ, IDUMMY,
485     &                 ISYOFZ, NLRZLBL, MAXZLBL, WORK, LWORK)
486      END IF
487
488C
489C------------------------------------------------------------------
490C     If solvent calculation then check vectors are solved
491C     self-consistently otherwise continue loop.
492C     CCSLV98, OC
493C------------------------------------------------------------------
494C
495      IF (CCSLV.OR.USE_PELIB()) THEN
496C
497         IF (CC3) CALL QUIT('No solvent for CC3...')
498C
499         LSLTCVG = .TRUE.
500         LSLLCVG = .TRUE.
501C
502C        Test if there exist vectors that are not converged.
503C
504         WRITE(LUPRI,'(//,1X,A)') 'Testing R1 list: '
505         WRITE(LUPRI,'(/,1X,A)')
506     * 'List nr.  Label     Symmetry      Frequency      T1-norm     '
507     * //'   T2-norm'
508         DO 997  ILSTNR = 1, NLRTLBL
509C
510            ISYM = ISYLRT(ILSTNR)
511            KT1  = 1
512            KT2  = KT1  + NT1AM(ISYM)
513            IF (.NOT. CCS) THEN
514               KEND1= KT2  + NT2AM(ISYM)
515            ELSE
516               KEND1= KT2
517            ENDIF
518            LEND1= LWORK - KEND1
519            IF (LEND1 .LE. 0 ) CALL QUIT(
520     *          'Insufficient workspace in CC_SOLRSP ')
521            IOPT  = 3
522            CALL CC_RDRSP('R1',ILSTNR,ISYM,IOPT,MODEL1,WORK(KT1),
523     *                    WORK(KT2))
524            X1  = DDOT(NT1AM(ISYM),WORK(KT1),1,WORK(KT1),1)
525            X2  = 0.0D0
526            IF (.NOT.CCS) THEN
527               X2  = DDOT(NT2AM(ISYM),WORK(KT2),1,WORK(KT2),1)
528            ENDIF
529            WRITE(LUPRI,'(1X,I5,5X,A8,2X,I4,6X,F15.10,F15.10,F15.10)')
530     *           ILSTNR,LRTLBL(ILSTNR),ISYM,FRQLRT(ILSTNR),X1,X2
531C
532            XLNCCCU = X1+X2
533            XLNCCPR = XLRT(ILSTNR)
534            IF (ABS(XLNCCPR-XLNCCCU).GT.CVGTSOL) LSLTCVG = .FALSE.
535            IF (IPRINT.GT.16)  WRITE(LUPRI,*) ' LSLTCVG = ',LSLTCVG
536            XLRT(ILSTNR) = XLNCCCU
537C
538 997     CONTINUE
539C
540         WRITE(LUPRI,'(//,1X,A)') 'Testing L1 list: '
541         WRITE(LUPRI,'(/,1X,A)')
542     * 'List nr.  Label     Symmetry      Frequency      T1-norm     '
543     * //'   T2-norm'
544         DO 998  ILSTNR = 1, NLRZLBL
545C
546            ISYM = ISYLRZ(ILSTNR)
547            KT1  = 1
548            KT2  = KT1  + NT1AM(ISYM)
549            IF (.NOT. CCS) THEN
550               KEND1= KT2  + NT2AM(ISYM)
551            ELSE
552               KEND1= KT2
553            ENDIF
554            LEND1= LWORK - KEND1
555            IF (LEND1 .LE. 0 ) CALL QUIT(
556     *          'Insufficient workspace in CC_SOLRSP ')
557            IOPT  = 3
558            CALL CC_RDRSP('L1',ILSTNR,ISYM,IOPT,MODEL1,WORK(KT1),
559     *                    WORK(KT2))
560            X1  = DDOT(NT1AM(ISYM),WORK(KT1),1,WORK(KT1),1)
561            X2  = 0.0D0
562            IF (.NOT.CCS) THEN
563               X2  = DDOT(NT2AM(ISYM),WORK(KT2),1,WORK(KT2),1)
564            ENDIF
565            WRITE(LUPRI,'(1X,I5,5X,A8,2X,I4,6X,F15.10,F15.10,F15.10)')
566     *           ILSTNR,LRZLBL(ILSTNR),ISYM,FRQLRZ(ILSTNR),X1,X2
567C
568            XLNCCCU = X1+X2
569            XLNCCPR = XLRZ(ILSTNR)
570            IF (ABS(XLNCCPR-XLNCCCU).GT.CVGLSOL) LSLLCVG = .FALSE.
571            IF (IPRINT.GT.16)  WRITE(LUPRI,*) ' LSLLCVG = ',LSLLCVG
572            XLRZ(ILSTNR) = XLNCCCU
573C
574 998     CONTINUE
575C
576C
577         IF (LSLTCVG.AND.LSLLCVG) THEN
578           WRITE(LUPRI,*)
579     *      ' Solvent 1-order response equations converged in'
580     *      ,ICCSLIT,' solvent iterations'
581         ELSE
582           ICCSLIT = ICCSLIT + 1
583           IF (ICCSLIT.GT.MXCCSLIT) THEN
584             CALL QUIT('Too many solvent iterations in CC_SOLRSP')
585           ENDIF
586! BUG FIX OVE DEC 2014
587           IF (CCTPA.AND.(.NOT.TPOLDW)) THEN
588              WRITE(LUPRI,*) 'NO SOLVENT ITERATION FOR TPA'
589           ELSE
590              WRITE(LUPRI,*)' Solvent 1-order response equations not ',
591     *                      'converged in',ICCSLIT,' solvent iterations'
592              WRITE(LUPRI,*)' A New iteration starts '
593              GOTO 999
594           END IF
595         ENDIF
596C
597      ENDIF
598C======================================================================
599C     Solve projected 1st-order Lagrange multiplier response equations.
600C======================================================================
601      LIST  = 'PL1'
602      NSTAT = 0
603      ORDER = 1
604      ISIDE = -1
605
606      IF ( .NOT. L1SKIP ) THEN
607        CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
608     &                 ISYSPL1, ISTPL1, EIGPL1, LPRPL1,
609     &                 ISYPL1,  LBLPL1, FRQPL1, IDUMMY,
610     &                 ISYOFPL1, NPL1LBL, MAXPL1LBL,
611     &                 WORK, LWORK)
612
613      ELSE
614        WRITE (LUPRI,*)
615     &        'Skipped proj. 1-order left response (PL1) equations'
616        L1SKIP = .FALSE.
617      END IF
618
619      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
620        CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
621     &                 ISYSPL1, ISTPL1, EIGPL1, LPRPL1,
622     &                 ISYPL1,  LBLPL1, FRQPL1, IDUMMY,
623     &                 ISYOFPL1, NPL1LBL, MAXPL1LBL, WORK, LWORK)
624      END IF
625
626C======================================================================
627C     Calculate F-transformed eigenvectors
628C======================================================================
629      LIST = 'FR '
630      NSTAT = 1
631      ORDER = 0
632      ISIDE = +1
633
634
635      IF ((.NOT.M1SKIP).AND.(.NOT.FRSKIP)) THEN
636
637         IF (NLRM .GT. 0) THEN
638
639            IF (NLRM .GT. MAXFDIM) THEN
640              WRITE (LUPRI,*)
641     &              'ERROR: Dimension MAXFDIM to small in CC_SOLRSP.'
642              CALL QUIT('Dimension MAXFDIM to small in CC_SOLRSP.')
643            END IF
644C
645C           set up list of F matrix transformations:
646C
647            DO IDXM = 1, NLRM
648              IDXRE = ILRM(IDXM)
649              IFTRAN(1,IDXM) = 0
650              IFTRAN(2,IDXM) = IDXRE
651              IFTRAN(3,IDXM) = IDXM
652            END DO
653C
654C           call CC_FMATRIX to do the F * RE transformations:
655C
656            IOPT = 3
657            CALL CC_FMATRIX(IFTRAN, NLRM, 'L0 ', 'RE ', IOPT, 'FR ',
658     &                      IDUMMY, DUMMY, 0, WORK, LWORK )
659
660
661         END IF
662
663      ELSE
664        WRITE (LUPRI,*) 'Skipped of F x RE transformations'
665        FRSKIP = .FALSE.
666      END IF
667
668      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
669        CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
670     &                 ISYLRM, ILRM, FRQLRM, LDUMMY,
671     &                 IDUMMY,   IDUMMY, DUMMY, IDUMMY,
672     &                 ISYOFM, NLRM, MAXM, WORK, LWORK)
673      END IF
674
675C======================================================================
676C     Solve equations for excited state M Lagrange multipliers:
677C======================================================================
678      LIST = 'M1 '
679      NSTAT = 1
680      ORDER = 0
681      ISIDE = -1
682
683      IF ( .NOT. M1SKIP ) THEN
684        CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
685     &                 ISYLRM, ILRM, FRQLRM, LDUMMY,
686     &                 IDUMMY,   IDUMMY, DUMMY, IDUMMY,
687     &                 ISYOFM, NLRM, MAXM,
688     &                 WORK, LWORK)
689      ELSE
690        WRITE (LUPRI,*) 'Skipped first-order M vector (M1) equations'
691        M1SKIP = .FALSE.
692      END IF
693
694      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
695        CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
696     &                 ISYLRM, ILRM, FRQLRM, LDUMMY,
697     &                 IDUMMY,   IDUMMY, DUMMY, IDUMMY,
698     &                 ISYOFM, NLRM, MAXM, WORK, LWORK)
699      END IF
700
701C======================================================================
702C     Make rhs for excited state Lagrangian multipliers E0:
703C======================================================================
704
705      LIST = 'BE '
706      NSTAT = 0
707      ORDER = 0
708      ISIDE = -1
709
710      IF (.NOT.BESKIP) THEN
711
712         IF (NXGRST .GT. 0) THEN
713
714            IF (NXGRST .GT. MAXFDIM)
715     &        CALL QUIT('Dimension MAXFDIM to small in CC_SOLRSP.')
716
717            ! set up list of 'F matrix' transformations:
718            DO IDXF = 1, NXGRST
719              IFTRAN(1,IDXF) = IXGRST(IDXF)
720              IFTRAN(2,IDXF) = IXGRST(IDXF)
721              IFTRAN(3,IDXF) = IDXF
722            END DO
723
724            ! call CC_FMATRIX to do the LE * B * RE transformations:
725            IOPT = 3
726            CALL CC_FMATRIX(IFTRAN, NXGRST, 'LE ', 'RE ', IOPT, LIST,
727     &                      IDUMMY, DUMMY, 0, WORK, LWORK )
728
729            ! add eta^(0) to complete rhs vectors:
730            NTAMP = NT1AMX
731            IF (.NOT.CCS) NTAMP = NTAMP + NT2AMX
732            IF (CCR12) NTAMP = NTAMP + NTR12AM(1)
733            KETA  = 1
734            KEND1 = KETA  + NTAMP
735            LWRK1 = LWORK - KEND1
736
737            KBEVC = KEND1
738            KEND2 = KBEVC + NTAMP
739            LWRK2 = LWORK - KEND2
740            IF (LWRK2.LT.0) CALL QUIT('Out of memory in CC_SOLRSP.')
741
742            CALL CC_ETA(WORK(KETA),WORK(KEND1),LWRK1)
743
744            DO IDXF = 1, NXGRST
745              IOPT = 3
746              CALL CC_RDRSP('BE ',IDXF,ISYM0,IOPT,MODEL,
747     &                       WORK(KBEVC),WORK(KBEVC+NT1AMX))
748              CALL DAXPY(NTAMP,ONE,WORK(KETA),1,WORK(KBEVC),1)
749              IOPT = 3
750              CALL CC_WRRSP('BE ',IDXF,ISYM0,IOPT,MODEL,DUMMY,
751     &               WORK(KBEVC),WORK(KBEVC+NT1AMX),WORK(KEND2),LWRK2)
752              IF (CC3) THEN
753                ! for CC3 also correct the effective RHS vectors
754                IOPT = 24
755                CALL CC_RDRSP('BE ',IDXF,ISYM0,IOPT,MODEL,
756     &                        WORK(KBEVC),WORK(KBEVC+NT1AMX))
757                CALL DAXPY(NTAMP,ONE,WORK(KETA),1,WORK(KBEVC),1)
758                IOPT = 24
759                CALL CC_WRRSP('BE ',IDXF,ISYM0,IOPT,MODEL,DUMMY,
760     &                WORK(KBEVC),WORK(KBEVC+NT1AMX),WORK(KEND2),LWRK2)
761              END IF
762            END DO
763         END IF
764
765      ELSE
766        WRITE (LUPRI,*) 'Skipped of LE x B x RE transformations'
767        BESKIP = .FALSE.
768      END IF
769
770      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
771        NSTAT = 1
772        CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
773     &                 IDUMMY,  IXGRST, EIGVAL, LDUMMY,
774     &                 IDUMMY,  IDUMMY, DUMMY, IDUMMY,
775     &                 ISYOFXG, NXGRST, MXXGST, WORK, LWORK)
776      END IF
777
778C======================================================================
779C     Solve equations for excited state lagrange multipliers E0:
780C======================================================================
781
782      LIST = 'E0 '
783      NSTAT = 0
784      ORDER = 0
785      ISIDE = -1
786
787      IF ( .NOT. E0SKIP ) THEN
788        CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
789     &                 IDUMMY,  IDUMMY, DUMMY,  LDUMMY,
790     &                 IDUMMY,  IDUMMY, DUMMY,  IDUMMY,
791     &                 ISYOFXG, NXGRST, MXXGST,
792     &                 WORK, LWORK)
793      ELSE
794        WRITE (LUPRI,*) 'Skipped exc. st. Lagr. mult.  (E0) equations'
795        E0SKIP = .FALSE.
796      END IF
797
798      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
799        NSTAT = 1
800        CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
801     &                 IDUMMY,  IXGRST, EIGVAL, LDUMMY,
802     &                 IDUMMY,  IDUMMY, DUMMY,  IDUMMY,
803     &                 ISYOFXG, NXGRST, MXXGST, WORK, LWORK)
804      END IF
805
806C======================================================================
807C     Make rhs for N2-vector equations.
808C======================================================================
809
810      LIST = 'BR '
811      NSTAT = 2
812      ORDER = 0
813      ISIDE = -1
814
815      IF ((.NOT.N2SKIP).AND.(.NOT.BRSKIP)) THEN
816
817         IF (NQRN2 .GT. 0) THEN
818
819            IF (NQRN2 .GT. MAXFDIM) THEN
820              WRITE (LUPRI,*)
821     &              'ERROR: Dimension MAXFDIM to small in CC_SOLRSP.'
822              CALL QUIT('Dimension MAXFDIM to small in CC_SOLRSP.')
823            END IF
824C
825C           set up list of 'F matrix' transformations:
826C
827            DO IDXF = 1, NQRN2
828              IFTRAN(1,IDXF) = IIN2(IDXF)
829              IFTRAN(2,IDXF) = IFN2(IDXF)
830              IFTRAN(3,IDXF) = IDXF
831            END DO
832C
833C           call CC_FMATRIX to do the LE * B * RE transformations:
834C
835            IOPT = 3
836            CALL CC_FMATRIX(IFTRAN, NQRN2, 'LE ', 'RE ', IOPT, LIST,
837     &                      IDUMMY, DUMMY, 0, WORK, LWORK )
838
839         END IF
840
841      ELSE
842        WRITE (LUPRI,*) 'Skipped of LE x B x RE transformations'
843        BRSKIP = .FALSE.
844      END IF
845
846      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
847        CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
848     &                 ISYSN2, ISTN2, EIGN2, LDUMMY,
849     &                 IDUMMY,   IDUMMY, DUMMY, IDUMMY,
850     &                 ISYOFN2, NQRN2, MAXQRN2, WORK, LWORK)
851      END IF
852
853C======================================================================
854C     Solve N2-vector for quadratic response 2. residue calculation.
855C======================================================================
856
857      LIST = 'N2 '
858      NSTAT = 2
859      ORDER = 0
860      ISIDE = -1
861
862      IF ( .NOT. N2SKIP ) THEN
863        CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
864     &                 ISYSN2, ISTN2, EIGN2, LDUMMY,
865     &                 IDUMMY,   IDUMMY, DUMMY, IDUMMY,
866     &                 ISYOFN2, NQRN2, MAXQRN2,
867     &                 WORK, LWORK)
868      ELSE
869        WRITE (LUPRI,*) 'Skipped second-order N vector (N2) equations'
870        N2SKIP = .FALSE.
871      END IF
872
873      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
874        CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
875     &                 ISYSN2, ISTN2, EIGN2, LDUMMY,
876     &                 IDUMMY,   IDUMMY, DUMMY, IDUMMY,
877     &                 ISYOFN2, NQRN2, MAXQRN2, WORK, LWORK)
878      END IF
879
880C======================================================================
881C     Solve first-order right Cauchy equations.
882C======================================================================
883      LIST  = 'RC '
884      NSTAT = 0
885      ORDER = 1
886      ISIDE = +1
887
888      IF ( .NOT. RCSKIP ) THEN
889
890       IF (NEWCAU) THEN
891
892        ! make sure that RC files are initialized:
893        ! if no restart vector available on file,
894        ! use next lower order RC vect.
895        DO IDXRC = 1, NLRCLBL
896          ISYM  = ISYLRC(IDXRC)
897          KVEC1 = 1
898          KVEC2 = KVEC1 + NT1AM(ISYM)
899          KEND1 = KVEC2 + NT2AM(ISYM)
900          LEND1 = LWORK - KEND1
901          IF (LEND1 .LT. 0) THEN
902            CALL QUIT('INSUFFICIENT WORKSPACE IN CC_SOLRSP.')
903          END IF
904          IOPT = 33
905          CALL CC_RDRSP('RC',IDXRC,ISYM,IOPT,MODEL,
906     &                  WORK(KVEC1),WORK(KVEC2)  )
907
908          IF (IOPT.EQ.33) THEN
909
910            NCAU = ILRCAU(IDXRC)
911            IF (NCAU.EQ.1) THEN
912              IDXR1 = IR1TAMP(LRCLBL(IDXRC),.FALSE.,0.0d0,ISYM)
913              CALL CC_RDRSP('R1',IDXR1,ISYM,IOPT,MODEL,
914     &                      WORK(KVEC1),WORK(KVEC2)  )
915            ELSE
916              IDXR1 = ILRCAMP(LRCLBL(IDXRC),NCAU-1,ISYM)
917              CALL CC_RDRSP('RC',IDXR1,ISYM,IOPT,MODEL,
918     &                      WORK(KVEC1),WORK(KVEC2)  )
919            END IF
920
921            CALL CC_WRRSP('RC',IDXRC,ISYM,IOPT,MODEL,DUMMY,
922     &                    WORK(KVEC1),WORK(KVEC2),WORK(KEND1),LEND1)
923          END IF
924
925        END DO
926       END IF
927
928       ! Put zeroth-order Cauchy moments at the end of the list,
929       ! but without increasing the number of equations that will
930       ! solved. The read routine CC_RDRSP will matp these vectors
931       ! onto the first-order amplitude response vectors which are
932       ! available from the R1 section above.
933       NLRCLBLSAVE = NLRCLBL
934       LRC1OPN = .TRUE. ! open list
935       DO IVEC = 1, NLRCLBLSAVE
936         INUM = ILRCAMP(LRCLBL(IVEC),0,ISYLRC(IVEC))
937       END DO
938       LRC1OPN = .FALSE. ! close list
939
940
941       CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
942     &                IDUMMY, IDUMMY, DUMMY, LDUMMY,
943     &                ISYLRC, LRCLBL, DUMMY, ILRCAU,
944     &                ISYOFC, NLRCLBLSAVE, MAXCLBL,
945     &                WORK, LWORK)
946
947      ELSE
948        WRITE (LUPRI,*) 'Skipped first-order right Cauchy equations'
949        RCSKIP = .FALSE.
950      END IF
951
952C======================================================================
953C     Add (unrelaxed) static R1 vectors to Cauchy RC list as RC(0)
954C======================================================================
955      LIST  = 'RC '
956      NSTAT = 0
957      ORDER = 1
958      ISIDE = +1
959
960      LRC1OPN = .TRUE.
961      DO IVEC = 1, NLRTLBL
962        IF (ABS(FRQLRT(IVEC)).LT.TOL .AND. (.NOT.LORXLRT(IVEC)) )  THEN
963          INUM = ILRCAMP(LRTLBL(IVEC),0,ISYLRT(IVEC))
964        END IF
965      END DO
966      LRC1OPN = .FALSE.
967
968      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
969       CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
970     &                IDUMMY, IDUMMY, DUMMY, LDUMMY,
971     &                ISYLRC, LRCLBL, DUMMY, ILRCAU,
972     &                ISYOFC, NLRCLBL, MAXCLBL, WORK, LWORK)
973      END IF
974
975C======================================================================
976C     precalculate F * RC transformations required to solve for
977C     the first-order left Cauchy vectors LC
978C======================================================================
979      LIST  = 'FC '
980      NSTAT = 0
981      ORDER = 1
982      ISIDE = +1
983
984      IF ( .NOT. FCSKIP ) THEN
985
986         IF (NLC1LBL .GT. 0) THEN
987C
988C          check compatibility with dimension of IFTRAN:
989C
990           IF (NLC1LBL .GT. MAXFDIM) THEN
991             WRITE (LUPRI,*)
992     &             'ERROR: Dimension MAXFDIM to small in CC_SOLRSP.'
993             CALL QUIT('Dimension MAXFDIM to small in CC_SOLRSP.')
994           END IF
995C
996C          set up list of F matrix transformations:
997C
998           DO IDXLC = 1, NLC1LBL
999             IDXRC=ILRCAMP(LBLLC1(IDXLC),ILC1CAU(IDXLC),ISYLC1(IDXLC))
1000             IFTRAN(1,IDXLC) = 0
1001             IFTRAN(2,IDXLC) = IDXRC
1002             IFTRAN(3,IDXLC) = IDXRC
1003           END DO
1004C
1005C          call CC_FMATRIX to do the F * RC transformations:
1006C
1007           IOPT = 3
1008           CALL CC_FMATRIX(IFTRAN, NLC1LBL, 'L0 ', 'RC ', IOPT, 'FC ',
1009     &                     IDUMMY, DUMMY, 0, WORK, LWORK )
1010
1011         ENDIF
1012
1013      ELSE
1014        WRITE (LUPRI,*) 'Skipped calculation of F x RC transformations'
1015        FCSKIP = .FALSE.
1016      END IF
1017
1018C==================================================================
1019C     Solve first-order left Cauchy equations:
1020C==================================================================
1021      LIST  = 'LC '
1022      NSTAT = 0
1023      ORDER = 1
1024      ISIDE = -1
1025
1026      IF ( .NOT. LCSKIP ) THEN
1027        CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
1028     &                 IDUMMY, IDUMMY, DUMMY, LDUMMY,
1029     &                 ISYLC1, LBLLC1, DUMMY, ILC1CAU,
1030     &                 ISYOFLC1, NLC1LBL, MAXLC1LBL,
1031     &                 WORK, LWORK)
1032      ELSE
1033        WRITE (LUPRI,*) 'Skipped first-order left Cauchy equations'
1034        LCSKIP = .FALSE.
1035      END IF
1036
1037C======================================================================
1038C     Add static L1 vectors to Cauchy LC list as LC(0)
1039C======================================================================
1040      LIST  = 'LC '
1041      NSTAT = 0
1042      ORDER = 1
1043      ISIDE = -1
1044
1045      LLC1OPN = .TRUE.
1046      DO IVEC = 1, NLRZLBL
1047        IF (ABS(FRQLRZ(IVEC)).LT.TOL .AND. (.NOT.LORXLRZ(IVEC)) )  THEN
1048          INUM = ILC1AMP(LRZLBL(IVEC),0,ISYLRZ(IVEC))
1049        END IF
1050      END DO
1051      LLC1OPN = .FALSE.
1052
1053      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
1054        CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
1055     &                 IDUMMY, IDUMMY, DUMMY, LDUMMY,
1056     &                 ISYLC1, LBLLC1, DUMMY, ILC1CAU,
1057     &                 ISYOFLC1, NLC1LBL, MAXLC1LBL, WORK, LWORK)
1058      END IF
1059
1060C======================================================================
1061C     precalculate right hand side vectors for second-order
1062C     cluster amplitude response equations.
1063C======================================================================
1064      LIST  = 'O2 '
1065      NSTAT = 0
1066      ORDER = 2
1067      ISIDE = +1
1068
1069      IF ( .NOT. O2SKIP ) THEN
1070
1071      CALL CCRHSVEC('O2 ',LBLO2,IDUMMY,IDUMMY,DUMMY,ISYO2,FRQO2,
1072     &              LORXO2,IDUMMY,NO2LBL,MAXO2LBL,0,
1073     &              WORK,LWORK)
1074
1075      ELSE
1076        WRITE (LUPRI,*)
1077     &        'Skipped calculation of O2 right hand side vectors.'
1078        O2SKIP = .FALSE.
1079      END IF
1080
1081      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
1082        CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
1083     &                 IDUMMY, IDUMMY, DUMMY, LDUMMY,
1084     &                 ISYO2, LBLO2, FRQO2, IDUMMY,
1085     &                 ISYOFO2, NO2LBL, MAXO2LBL, WORK, LWORK)
1086      END IF
1087
1088C======================================================================
1089C     Solve second-order right cluster amplitude response equations.
1090C======================================================================
1091      LIST  = 'R2 '
1092      NSTAT = 0
1093      ORDER = 2
1094      ISIDE = +1
1095
1096      IF ( .NOT. R2SKIP ) THEN
1097        CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
1098     &                 IDUMMY, IDUMMY, DUMMY, LDUMMY,
1099     &                 ISYR2T, LBLR2T, FRQR2T, IDUMMY,
1100     &                 ISYOFT2, NR2TLBL, MAXT2LBL,
1101     &                 WORK, LWORK)
1102      ELSE
1103        WRITE (LUPRI,*) 'Skipped second-order response (R2) equations'
1104        R2SKIP = .FALSE.
1105      END IF
1106
1107      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
1108        CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
1109     &                 IDUMMY, IDUMMY, DUMMY, LDUMMY,
1110     &                 ISYR2T, LBLR2T, FRQR2T, IDUMMY,
1111     &                 ISYOFT2, NR2TLBL, MAXT2LBL, WORK, LWORK)
1112      END IF
1113
1114C======================================================================
1115C     precalculate CO2 vectors which contain all contributions to the
1116C     rhs vectors for the second-order right Cauchy vectors which are
1117C     independent of second-order right Cauchy vectors
1118C======================================================================
1119      LIST  = 'CO2'
1120      NSTAT = 0
1121      ORDER = 2
1122      ISIDE = +1
1123
1124      IF ( .NOT. CO2SKIP ) THEN
1125
1126        CALL CCRHSVEC('CO2',LBLCO2,IDUMMY,IDUMMY,DUMMY,ISYCO2,DUMMY,
1127     &                LDUMMY,ICO2CAU,NCO2LBL,MAXCO2LBL,0,
1128     &                WORK,LWORK)
1129
1130      ELSE
1131        WRITE (LUPRI,*)
1132     &        'Skipped calculation of CO2 right hand side vectors.'
1133        CO2SKIP = .FALSE.
1134      END IF
1135
1136C===================================================================
1137C     Solve second-order right Cauchy response equations.
1138C===================================================================
1139      LIST  = 'CR2'
1140      NSTAT = 0
1141      ORDER = 2
1142      ISIDE = +1
1143
1144      IF ( .NOT. CR2SKIP ) THEN
1145        CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
1146     &                 IDUMMY, IDUMMY, DUMMY, LDUMMY,
1147     &                 ISYCR2, LBLCR2, DUMMY, ICR2CAU,
1148     &                 ISYOFCR2, NCR2LBL, MAXCR2LBL,
1149     &                 WORK, LWORK)
1150      ELSE
1151        WRITE (LUPRI,*)
1152     &        'Skipped second-order right Cauchy (CR2) equations'
1153        CR2SKIP = .FALSE.
1154      END IF
1155
1156C======================================================================
1157C     Add static R2 vectors to second-order Cauchy CR2 list as CR2(0,0)
1158C======================================================================
1159      LIST  = 'CR2'
1160      NSTAT = 0
1161      ORDER = 2
1162      ISIDE = +1
1163
1164      LCR2OPN = .TRUE.
1165      DO IVEC = 1, NR2TLBL
1166        IF (ABS(FRQR2T(IVEC,1)).LT.TOL .AND. ABS(FRQR2T(IVEC,2)).LT.TOL
1167     &    .AND. (.NOT.LORXR2T(IVEC,1)) .AND. (.NOT.LORXR2T(IVEC,2))
1168     &     )  THEN
1169          INUM = ICR2AMP(LBLR2T(IVEC,1),0,ISYR2T(IVEC,1),
1170     &                   LBLR2T(IVEC,2),0,ISYR2T(IVEC,2))
1171        END IF
1172      END DO
1173      LCR2OPN = .FALSE.
1174
1175      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
1176        CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
1177     &                 IDUMMY, IDUMMY, DUMMY, LDUMMY,
1178     &                 ISYCR2, LBLCR2, DUMMY, ICR2CAU,
1179     &                 ISYOFCR2, NCR2LBL, MAXCR2LBL, WORK, LWORK)
1180      END IF
1181
1182C======================================================================
1183C     Add static O2 vectors to second-order Cauchy CO2 list as CO2(0,0)
1184C======================================================================
1185      LIST  = 'CO2'
1186      NSTAT = 0
1187      ORDER = 2
1188      ISIDE = +1
1189
1190      LCO2OPN = .TRUE.
1191      DO IVEC = 1, NO2LBL
1192        IF (ABS(FRQO2(IVEC,1)).LT.TOL .AND. ABS(FRQO2(IVEC,2)).LT.TOL
1193     &     )  THEN
1194          INUM = IRHSCR2(LBLO2(IVEC,1),0,ISYO2(IVEC,1),
1195     &                   LBLO2(IVEC,2),0,ISYO2(IVEC,2))
1196        END IF
1197      END DO
1198      LCO2OPN = .FALSE.
1199
1200      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
1201        CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
1202     &                 IDUMMY, IDUMMY, DUMMY, LDUMMY,
1203     &                 ISYCO2, LBLCO2, DUMMY, ICO2CAU,
1204     &                 ISYOFCO2, NCO2LBL, MAXCO2LBL, WORK, LWORK)
1205      END IF
1206
1207C======================================================================
1208C     precalculate second-order chi vectors required to set up the
1209C     right hande side vector for the second-order multiplier equations
1210C======================================================================
1211      IF ( .NOT. X2SKIP ) THEN
1212
1213        CALL CC_ETADRV('X2 ',LBLX2,IDUMMY,IDUMMY,DUMMY,
1214     &                 ISYX2,FRQX2,IDUMMY,NX2LBL,MAXX2LBL,
1215     &                 WORK,LWORK)
1216
1217      ELSE
1218        WRITE (LUPRI,*)
1219     &        'Skipped calculation of X2 right hand side vectors.'
1220        X2SKIP = .FALSE.
1221      END IF
1222
1223C======================================================================
1224C     precalculate F * R2 transformations required to solve for
1225C     the second order coupled cluster multiplier reponses L2
1226C======================================================================
1227      IF ( .NOT. F2SKIP ) THEN
1228
1229         IF (NL2LBL .GT. 0) THEN
1230C
1231C          check compatibility with dimension of IFTRAN:
1232C
1233           IF (NL2LBL .GT. MAXFDIM) THEN
1234             WRITE (LUPRI,*)
1235     &             'ERROR: Dimension MAXFDIM to small in CC_SOLRSP.'
1236             CALL QUIT('Dimension MAXFDIM to small in CC_SOLRSP.')
1237           END IF
1238C
1239C          set up list of F matrix transformations:
1240C
1241           DO IDXL2 = 1, NL2LBL
1242            IDXR2 =
1243     &        IR2TAMP(LBLAL2(IDXL2),.FALSE.,FRQAL2(IDXL2),ISYAL2(IDXL2),
1244     &                LBLBL2(IDXL2),.FALSE.,FRQBL2(IDXL2),ISYBL2(IDXL2))
1245
1246            IFTRAN(1,IDXL2) = 0
1247            IFTRAN(2,IDXL2) = IDXR2
1248            IFTRAN(3,IDXL2) = IDXR2
1249           END DO
1250C
1251C          call CC_FMATRIX to do the F * R2 transformations:
1252C
1253           IOPT = 3
1254           CALL CC_FMATRIX(IFTRAN, NL2LBL, 'L0 ', 'R2 ', IOPT, 'F2 ',
1255     &                     IDUMMY, DUMMY, 0, WORK, LWORK )
1256
1257         ENDIF
1258
1259      ELSE
1260        WRITE (LUPRI,*) 'Skipped calculation of F x R2 transformations'
1261        F2SKIP = .FALSE.
1262      END IF
1263
1264C==================================================================
1265C     Solve second-order lagrange multiplier response equations.
1266C==================================================================
1267      LIST  = 'L2 '
1268      NSTAT = 0
1269      ORDER = 2
1270      ISIDE = -1
1271
1272      IF ( .NOT. L2SKIP ) THEN
1273        CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
1274     &                 IDUMMY, IDUMMY, DUMMY, LDUMMY,
1275     &                 ISYL2, LBLL2, FRQL2, IDUMMY,
1276     &                 ISYOFL2, NL2LBL, MAXL2LBL,
1277     &                 WORK, LWORK)
1278      ELSE
1279        WRITE (LUPRI,*)
1280     &        'Skipped second-order left response (L2) equations'
1281        L2SKIP = .FALSE.
1282      END IF
1283
1284      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
1285        CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
1286     &                 IDUMMY, IDUMMY, DUMMY, LDUMMY,
1287     &                 ISYL2, LBLL2, FRQL2, IDUMMY,
1288     &                 ISYOFL2, NL2LBL, MAXL2LBL, WORK, LWORK)
1289      END IF
1290
1291C======================================================================
1292C     precalculate second-order Cauchy eta vectors required to set up
1293C     the right hande side vectors for the second-order left Cauchy
1294C     equations
1295C======================================================================
1296      IF ( .NOT. CX2SKIP ) THEN
1297
1298      CALL CC_ETADRV('CX2',LBLCX2,IDUMMY,IDUMMY,DUMMY,
1299     &               ISYCX2,DUMMY,ICX2CAU,NCX2LBL,MAXCX2LBL,
1300     &               WORK,LWORK)
1301
1302      ELSE
1303        WRITE (LUPRI,*)
1304     &        'Skipped calculation of CX2 right hand side vectors.'
1305        CX2SKIP = .FALSE.
1306      END IF
1307C======================================================================
1308C     precalculate F * CR2 transformations required to solve for
1309C     the second order left Cauchy vectors CL2
1310C======================================================================
1311      IF ( .NOT. CF2SKIP ) THEN
1312
1313         IF (NCL2LBL .GT. 0) THEN
1314C
1315C          check compatibility with dimension of IFTRAN:
1316C
1317           IF (NCL2LBL .GT. MAXFDIM) THEN
1318             WRITE (LUPRI,*)
1319     &             'ERROR: Dimension MAXFDIM to small in CC_SOLRSP.'
1320             CALL QUIT('Dimension MAXFDIM to small in CC_SOLRSP.')
1321           END IF
1322C
1323C          set up list of F matrix transformations:
1324C
1325           DO IDXL2 = 1, NCL2LBL
1326             IDXR2 = ICR2AMP(
1327     &                LBLCL2(IDXL2,1),ICL2CAU(IDXL2,1),ISYCL2(IDXL2,1),
1328     &                LBLCL2(IDXL2,2),ICL2CAU(IDXL2,2),ISYCL2(IDXL2,2))
1329
1330             IFTRAN(1,IDXL2) = 0
1331             IFTRAN(2,IDXL2) = IDXR2
1332             IFTRAN(3,IDXL2) = IDXR2
1333           END DO
1334C
1335C          call CC_FMATRIX to do the F * CR2 transformations:
1336C
1337           IOPT = 3
1338           CALL CC_FMATRIX(IFTRAN, NCL2LBL, 'L0 ', 'CR2', IOPT, 'CF2',
1339     &                     IDUMMY, DUMMY, 0, WORK, LWORK )
1340
1341         ENDIF
1342
1343      ELSE
1344        WRITE (LUPRI,*) 'Skipped calculation of F x CR2 transformations'
1345        CF2SKIP = .FALSE.
1346      END IF
1347
1348C==================================================================
1349C     Solve second-order left Cauchy vector ("CL2") equations.
1350C==================================================================
1351      LIST  = 'CL2'
1352      NSTAT = 0
1353      ORDER = 2
1354      ISIDE = -1
1355
1356      IF ( .NOT. CL2SKIP ) THEN
1357        CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
1358     &                 IDUMMY, IDUMMY, DUMMY, LDUMMY,
1359     &                 ISYCL2, LBLCL2, DUMMY, ICL2CAU,
1360     &                 ISYOFCL2, NCL2LBL, MAXCL2LBL,
1361     &                 WORK, LWORK)
1362      ELSE
1363        WRITE (LUPRI,*)
1364     &        'Skipped second-order left Cauchy (CL2) equations'
1365        CL2SKIP = .FALSE.
1366      END IF
1367
1368C======================================================================
1369C     Add static L2 vectors to second-order Cauchy CL2 list as CL2(0,0)
1370C======================================================================
1371      LIST  = 'CL2'
1372      NSTAT = 0
1373      ORDER = 2
1374      ISIDE = -1
1375
1376      LCL2OPN = .TRUE.
1377      DO IVEC = 1, NL2LBL
1378        IF (ABS(FRQL2(IVEC,1)).LT.TOL .AND. ABS(FRQL2(IVEC,2)).LT.TOL
1379     &    .AND. (.NOT.LORXL2(IVEC,1)) .AND. (.NOT.LORXL2(IVEC,2))
1380     &     )  THEN
1381          INUM = ICL2AMP(LBLL2(IVEC,1),0,ISYL2(IVEC,1),
1382     &                   LBLL2(IVEC,2),0,ISYL2(IVEC,2))
1383        END IF
1384      END DO
1385      LCL2OPN = .FALSE.
1386
1387      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
1388        CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
1389     &                 IDUMMY, IDUMMY, DUMMY, LDUMMY,
1390     &                 ISYCL2, LBLCL2, DUMMY, ICL2CAU,
1391     &                 ISYOFCL2, NCL2LBL, MAXCL2LBL, WORK, LWORK)
1392      END IF
1393
1394C======================================================================
1395C     Add static X2 vectors to second-order Cauchy CX2 list as CX2(0,0)
1396C======================================================================
1397      LCX2OPN = .TRUE.
1398      DO IVEC = 1, NX2LBL
1399        IF (ABS(FRQX2(IVEC,1)).LT.TOL .AND. ABS(FRQX2(IVEC,2)).LT.TOL
1400     &     )  THEN
1401          INUM = IETACL2(LBLX2(IVEC,1),0,ISYX2(IVEC,1),
1402     &                   LBLX2(IVEC,2),0,ISYX2(IVEC,2))
1403        END IF
1404      END DO
1405      LCX2OPN = .FALSE.
1406
1407C===================================================================
1408C     precalculate right hand side vectors for third-order
1409C     cluster amplitude response equations.
1410C===================================================================
1411      CALL CCRHSVEC('O3 ',LBLO3,IDUMMY,IDUMMY,DUMMY,ISYO3,FRQO3,
1412     &              LORXO3,IDUMMY,NO3LBL,MAXO3LBL,0,
1413     &              WORK,LWORK)
1414
1415
1416C===================================================================
1417C     Solve third-order right cluster amplitude response equations.
1418C===================================================================
1419      LIST  = 'R3 '
1420      NSTAT = 0
1421      ORDER = 3
1422      ISIDE = +1
1423
1424      CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
1425     &               IDUMMY, IDUMMY, DUMMY, LDUMMY,
1426     &               ISYR3T, LBLR3T, FRQR3T, IDUMMY,
1427     &               ISYOFT3, NR3TLBL, MAXT3LBL,
1428     &               WORK, LWORK)
1429
1430
1431C======================================================================
1432C     precalculate third-order chi vectors required to set up the
1433C     right hande side vector for the third-order multiplier equations
1434C======================================================================
1435
1436      CALL CC_ETADRV('X3 ',LBLX3,IDUMMY,IDUMMY,DUMMY,
1437     &               ISYX3,FRQX3,IDUMMY,NX3LBL,MAXX3LBL,
1438     &               WORK,LWORK)
1439
1440C======================================================================
1441C     precalculate F * R3 transformations required to solve for
1442C     the second order coupled cluster multiplier reponses L3
1443C======================================================================
1444      IF (NL3LBL .GT. 0) THEN
1445C
1446C       check compatibility with dimension of IFTRAN:
1447C
1448        IF (NL3LBL .GT. MAXFDIM) THEN
1449          WRITE (LUPRI,*)
1450     &          'ERROR: Dimension MAXFDIM to small in CC_SOLRSP.'
1451          CALL QUIT('Dimension MAXFDIM to small in CC_SOLRSP.')
1452        END IF
1453C
1454C       set up list of F matrix transformations:
1455C
1456        DO IDXL3 = 1, NL3LBL
1457          IDXR3 = IR3TAMP(LBLL3(IDXL3,1),FRQL3(IDXL3,1),ISYL3(IDXL3,1),
1458     &                    LBLL3(IDXL3,2),FRQL3(IDXL3,2),ISYL3(IDXL3,2),
1459     &                    LBLL3(IDXL3,3),FRQL3(IDXL3,3),ISYL3(IDXL3,3))
1460
1461          IFTRAN(1,IDXL3) = 0
1462          IFTRAN(2,IDXL3) = IDXR3
1463          IFTRAN(3,IDXL3) = IDXR3
1464        END DO
1465C
1466C       call CC_FMATRIX to do the F * R3 transformations:
1467C
1468        IOPT = 3
1469        CALL CC_FMATRIX(IFTRAN, NL3LBL, 'L0 ', 'R3 ', IOPT, 'F3 ',
1470     &                  IDUMMY, DUMMY, 0, WORK, LWORK )
1471
1472      ENDIF
1473
1474C==================================================================
1475C     Solve second-order lagrange multiplier response equations.
1476C==================================================================
1477      LIST  = 'L3 '
1478      NSTAT = 0
1479      ORDER = 3
1480      ISIDE = -1
1481
1482      CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
1483     &               IDUMMY, IDUMMY, DUMMY, LDUMMY,
1484     &               ISYL3, LBLL3, FRQL3, IDUMMY,
1485     &               ISYOFL3, NL3LBL, MAXL3LBL,
1486     &               WORK, LWORK)
1487
1488
1489C===================================================================
1490C     precalculate right hand side vectors for first-order
1491C     right eigenvector response equations.
1492C===================================================================
1493      CALL CCRHSVEC('EO1',LBLER1,ISYSER1,ISTER1,EIGER1,ISYOER1,FRQER1,
1494     &              LORXER1,IDUMMY,NER1LBL,MAXER1LBL,0,
1495     &              WORK,LWORK)
1496
1497
1498C===================================================================
1499C     Solve first-order right eigenvector response equations.
1500C===================================================================
1501      LIST  = 'ER1'
1502      NSTAT = 1
1503      ORDER = 1
1504      ISIDE = +1
1505
1506      CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
1507     &               ISYSER1, ISTER1, EIGER1, LPRER1,
1508     &               ISYOER1, LBLER1, FRQER1, IDUMMY,
1509     &               ISYOFER1, NER1LBL, MAXER1LBL,
1510     &               WORK, LWORK)
1511
1512      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
1513        CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
1514     &                 ISYSER1, ISTER1, EIGER1, LPRER1,
1515     &                 ISYOER1, LBLER1, FRQER1, IDUMMY,
1516     &                 ISYOFER1, NER1LBL, MAXER1LBL, WORK, LWORK)
1517      END IF
1518
1519C===================================================================
1520C     precalculate right hand side vectors for second-order
1521C     right eigenvector response equations.
1522C===================================================================
1523      CALL CCRHSVEC('EO2',LBLER2,ISYSER2,ISTER2,EIGER2,ISYOER2,FRQER2,
1524     &              LORXER2,IDUMMY,NER2LBL,MAXER2LBL,0,
1525     &              WORK,LWORK)
1526
1527
1528C===================================================================
1529C     Solve second-order right eigenvector response equations.
1530C===================================================================
1531      LIST  = 'ER2'
1532      NSTAT = 1
1533      ORDER = 2
1534      ISIDE = +1
1535
1536      CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
1537     &               ISYSER2, ISTER2, EIGER2, LPRER2,
1538     &               ISYOER2, LBLER2, FRQER2, IDUMMY,
1539     &               ISYOFER2, NER2LBL, MAXER2LBL,
1540     &               WORK, LWORK)
1541
1542C===================================================================
1543C     precalculate right hand side vectors (EX1) for first-order
1544C     left eigenvector response (EL1) equations.
1545C===================================================================
1546      if (lold) then
1547      WRITE(LUPRI,*) 'CC_SOLRSP: use old CC_ETADRV for EX1...'
1548      CALL CC_ETADRV('EX1',LBLEL1,ISYSEL1,ISTEL1,EIGEL1,
1549     &             ISYOEL1,FRQEL1,IDUMMY,NEL1LBL,MAXEL1LBL,
1550     &             WORK,LWORK)
1551      else
1552       CALL CC_ETADRV1('EX1',LBLEL1,ISYSEL1,ISTEL1,EIGEL1,
1553     &               ISYOEL1,FRQEL1,LORXEL1,IDUMMY,NEL1LBL,MAXEL1LBL,
1554     &               WORK,LWORK)
1555      end if
1556
1557C===================================================================
1558C     Solve first-order left eigenvector response equations.
1559C===================================================================
1560      LIST  = 'EL1'
1561      NSTAT = 1
1562      ORDER = 1
1563      ISIDE = -1
1564
1565      CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
1566     &               ISYSEL1, ISTEL1, EIGEL1, LPREL1,
1567     &               ISYOEL1, LBLEL1, FRQEL1, IDUMMY,
1568     &               ISYOFEL1, NEL1LBL, MAXEL1LBL,
1569     &               WORK, LWORK)
1570
1571
1572      IF (DEBUG .OR. LOCDBG. OR. IPRINT.GE.5) THEN
1573        CALL CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
1574     &                 ISYSEL1, ISTEL1, EIGEL1, LPREL1,
1575     &                 ISYOEL1, LBLEL1, FRQEL1, IDUMMY,
1576     &                 ISYOFEL1, NEL1LBL, MAXEL1LBL, WORK, LWORK)
1577      END IF
1578
1579C===================================================================
1580C     precalculate right hand side vectors for second-order
1581C     left eigenvector response equations.
1582C===================================================================
1583      CALL CC_ETADRV('EX2',LBLEL2,ISYSEL2,ISTEL2,EIGEL2,
1584     &               ISYOEL2,FRQEL2,IDUMMY,NEL2LBL,MAXEL2LBL,
1585     &               WORK,LWORK)
1586
1587
1588C===================================================================
1589C     Solve first-order left eigenvector response equations.
1590C===================================================================
1591      LIST  = 'EL2'
1592      NSTAT = 1
1593      ORDER = 2
1594      ISIDE = -1
1595
1596      CALL CC_SOLDRV(LIST,NSTAT,ORDER,ISIDE,APROXR12,
1597     &               ISYSEL2, ISTEL2, EIGEL2, LPRER2,
1598     &               ISYOEL2, LBLEL2, FRQEL2, IDUMMY,
1599     &               ISYOFEL2, NEL2LBL, MAXEL2LBL,
1600     &               WORK, LWORK)
1601
1602C---------------------------------------------------------------------
1603C     TEST - Loop over list of response vectors
1604C---------------------------------------------------------------------
1605C
1606C     Test average values for freq = 0 and tot sym vectors.
1607C
1608      IF (DEBUG .OR. (IPRINT.GE.5)) THEN
1609         DO ILSTNR = 1, NLRTLBL
1610C
1611            ISYM = ISYLRT(ILSTNR)
1612            FREK = FRQLRT(ILSTNR)
1613          IF (ISYM.EQ.1 .AND. DABS(FREK).LE.1.0D-10) THEN
1614C
1615            KT1  = 1
1616            KT2  = KT1  + NT1AM(ISYM)
1617            IF (.NOT. CCS) THEN
1618               KEND1= KT2  + NT2AM(ISYM)
1619            ELSE
1620               KEND1= KT2
1621            ENDIF
1622            IF (CCR12) THEN
1623               KTR12= KEND1
1624               KEND1= KTR12 + NTR12AM(ISYM)
1625            ENDIF
1626            LEND1= LWORK - KEND1
1627            IF (LEND1 .LE. 0 ) CALL QUIT(
1628     *       'Insufficient workspace in cc_solrsp ')
1629            IOPT  = 3
1630            CALL CC_RDRSP('R1 ',ILSTNR,ISYM,IOPT,MODEL1,WORK(KT1),
1631     *                    WORK(KT2))
1632            IF (CCR12) THEN
1633              IOPT = 32
1634              CALL CC_RDRSP('R1 ',ILSTNR,ISYM,IOPT,MODEL1,DUMMY,
1635     *                    WORK(KTR12))
1636            ENDIF
1637            CALL CC_TSTAV(ILSTNR,WORK(KT1),WORK(KEND1),LEND1,0)
1638C
1639          END IF
1640
1641         END DO
1642      ENDIF
1643
1644      IF (DEBUG .OR. (IPRINT.GE.5)) THEN
1645         DO ILSTNR = 1, NO1LBL
1646C
1647            ISYM = ISYO1(ILSTNR)
1648            FREK = FRQO1(ILSTNR)
1649          IF (ISYM.EQ.1 .AND. DABS(FREK).LE.1.0D-10) THEN
1650C
1651            KT1  = 1
1652            KT2  = KT1  + NT1AM(ISYM)
1653            IF (.NOT. CCS) THEN
1654               KEND1= KT2  + NT2AM(ISYM)
1655            ELSE
1656               KEND1= KT2
1657            ENDIF
1658            IF (CCR12) THEN
1659               KTR12= KEND1
1660               KEND1= KTR12 + NTR12AM(ISYM)
1661            ENDIF
1662            LEND1= LWORK - KEND1
1663            IF (LEND1 .LE. 0 ) CALL QUIT(
1664     *       'Insufficient workspace in cc_solrsp ')
1665            IOPT  = 3
1666            IF (CCSDT) IOPT = 24
1667            CALL CC_RDRSP('O1 ',ILSTNR,ISYM,IOPT,MODEL1,WORK(KT1),
1668     *                    WORK(KT2))
1669            IF (CCR12) THEN
1670              IOPT = 32
1671              CALL CC_RDRSP('O1 ',ILSTNR,ISYM,IOPT,MODEL1,DUMMY,
1672     *                    WORK(KTR12))
1673            ENDIF
1674            CALL CC_TSTAV(ILSTNR,WORK(KT1),WORK(KEND1),LEND1,1)
1675C
1676          END IF
1677
1678         END DO
1679      ENDIF
1680C
1681      IF (DEBUG .OR. (IPRINT.GE.5)) THEN
1682         DO ILSTNR = 1, NR2TLBL
1683            ISYM  = ILSTSYM('R2 ',ILSTNR)
1684            KT1   = 1
1685            KT2   = KT1  + NT1AM(ISYM)
1686            KEND1 = KT2
1687            IF (.NOT. CCS) KEND1 = KT2  + NT2AM(ISYM)
1688            IF (CCR12) THEN
1689               KTR12= KEND1
1690               KEND1= KTR12 + NTR12AM(ISYM)
1691            ENDIF
1692            LEND1 = LWORK - KEND1
1693            IF (LEND1.LE.0) CALL QUIT('Out of workspace in cc_solrsp')
1694            IOPT  = 3
1695            CALL CC_RDRSP('R2 ',ILSTNR,ISYM,IOPT,MODEL1,WORK(KT1),
1696     *                    WORK(KT2))
1697            IF (CCR12) THEN
1698              IOPT = 32
1699              CALL CC_RDRSP('R2 ',ILSTNR,ISYM,IOPT,MODEL1,DUMMY,
1700     *                    WORK(KTR12))
1701            ENDIF
1702            CALL CC_TSTAV2(ILSTNR,WORK(KT1),WORK(KEND1),LEND1,0)
1703         END DO
1704      ENDIF
1705
1706C
1707      WRITE (LUPRI,*)
1708      WRITE (LUPRI,*) 'Solution of CC response equations completed.'
1709      WRITE (LUPRI,*)
1710      CALL FLSHFO(LUPRI)
1711C
1712      CALL QEXIT('CC_SOLRSP')
1713      RETURN
1714C
1715   2  FORMAT ((1X,A1,1X,A4,1X,A1,1X,I4,3X,A1,1X,A52,1X,A1))
1716      END
1717*=====================================================================*
1718      SUBROUTINE CC_TSTLST(LIST,NSTAT,ORDER,ISIDE,
1719     &                     ISYMS, ISTAT, EIGVAL, LPROJ,
1720     &                     ISYMO, LABEL, FREQS,  ICAU,
1721     &                     ISYOF, NVEC,  MAXVEC, WORK, LWORK)
1722*---------------------------------------------------------------------*
1723*
1724*     Purpose: calculate and print norm of all vectors on a list
1725*              (mainly for testing and debugging purposes)
1726*
1727*     Christof Haettig, August 1999
1728*
1729*---------------------------------------------------------------------*
1730      IMPLICIT NONE
1731C
1732#include "priunit.h"
1733#include "maxorb.h"
1734#include "ccorb.h"
1735#include "ccsdsym.h"
1736#include "ccsdinp.h"
1737#include "dummy.h"
1738Cholesky
1739#include "ccdeco.h"
1740      LOGICAL CCSEFF
1741Cholesky
1742C
1743      LOGICAL TRIPLET
1744      PARAMETER (TRIPLET= .FALSE.)
1745C
1746      CHARACTER LIST*(3)
1747
1748      INTEGER NSTAT, ORDER, MAXVEC, NVEC, ISIDE, LWORK
1749      INTEGER ISYOF(8), KT12, NCCVAR
1750      INTEGER ISYM, KT1, KT2, IOPT, ILSTNR, KEND1, LEND1, NTAMP
1751
1752      LOGICAL LPROJ(MAXVEC)
1753      CHARACTER*8 LABEL(MAXVEC,ORDER)
1754      CHARACTER*10 MODEL
1755      INTEGER ISYMS(MAXVEC,NSTAT), ISTAT(MAXVEC,NSTAT)
1756      INTEGER ISYMO(MAXVEC,ORDER), ICAU(MAXVEC,ORDER)
1757
1758#if defined (SYS_CRAY)
1759      REAL FREQS(MAXVEC,ORDER), EIGVAL(MAXVEC,NSTAT)
1760      REAL WORK(LWORK)
1761      REAL X1, X2, X12, DDOT
1762#else
1763      DOUBLE PRECISION FREQS(MAXVEC,ORDER), EIGVAL(MAXVEC,NSTAT)
1764      DOUBLE PRECISION WORK(LWORK)
1765      DOUBLE PRECISION X1, X2, X12, DDOT
1766#endif
1767
1768      INTEGER ILSTSYM
1769
1770      CCSEFF = CCS .OR. (CHOINT .AND. CC2)  ! No doubles for Chol. CC2, tbp
1771
1772
1773C
1774         IF (NVEC.LE.0) THEN
1775           WRITE (LUPRI,'(//1X,3A)')  'List ',LIST,' is empty.'
1776           RETURN
1777         END IF
1778C
1779         IF     (LIST(2:2).NE.'C'.AND.ORDER.EQ.1.AND.NSTAT.EQ.0) THEN
1780
1781           WRITE (LUPRI,'(//1X,3A,I3,A,/1X,A)') 'Testing ',LIST,
1782     &        ' list with',NVEC, ' elements on it:',
1783     &        'List Nr.  Label     Symmetry      Frequency'
1784     &      //"      T1-norm        T2-norm     T2'-norm"
1785
1786         ELSEIF (LIST(1:1).NE.'C'.AND.ORDER.EQ.2.AND.NSTAT.EQ.0) THEN
1787
1788           WRITE (LUPRI,'(//1X,3A,I3,A,/1X,A)') 'Testing ',LIST,
1789     &        ' list with', NVEC, ' elements on it:',
1790     &        'List Nr.  Label 1   Label 2   Symmetry      Freq. 1'
1791     &      //'     Freq. 2      T1-norm        T2-norm'
1792
1793         ELSEIF (LIST(2:2).EQ.'C'.AND.ORDER.EQ.1.AND.NSTAT.EQ.0) THEN
1794
1795           WRITE (LUPRI,'(//1X,3A,I3,A,/1X,A)') 'Testing ',LIST,
1796     &        ' list with', NVEC, ' elements on it:',
1797     &        'List Nr.  Label     Symmetry     Order'
1798     &      //'    T1-norm             T2-norm'
1799
1800         ELSEIF (LIST(1:1).EQ.'C'.AND.ORDER.EQ.2.AND.NSTAT.EQ.0) THEN
1801
1802           WRITE (LUPRI,'(//1X,3A,I3,A,/1X,A)') 'Testing ',LIST,
1803     &        ' list with', NVEC, ' elements on it:',
1804     &        'List Nr.  Label 1   Label 2   Symmetry   Order 1'
1805     &      //'   Order 2  T1-norm             T2-norm'
1806
1807         ELSEIF (ORDER.EQ.0 .AND. NSTAT.EQ.1) THEN
1808
1809           WRITE (LUPRI,'(//1X,3A,I3,A,/1X,A)') 'Testing ',LIST,
1810     &        ' list with', NVEC, ' elements on it:',
1811     &        'List Nr.  State     Symmetry     Eigenvalue'
1812     &      //'      T1-norm        T2-norm'
1813
1814         ELSEIF (ORDER.EQ.1 .AND. NSTAT.EQ.1) THEN
1815
1816           WRITE (LUPRI,'(//1X,3A,I3,A,/1X,A)') 'Testing ',LIST,
1817     &        ' list with', NVEC, ' elements on it:',
1818     &        'List Nr.  Label     State     Symmetry     Eigenval.'
1819     &      //'   Frequency    T1-norm        T2-norm'
1820
1821         ELSEIF (ORDER.EQ.0 .AND. NSTAT.EQ.2) THEN
1822
1823           WRITE (LUPRI,'(//1X,3A,I3,A,/1X,A)') 'Testing ',LIST,
1824     &        ' list with', NVEC, ' elements on it:',
1825     &        'List Nr.  State 1   State 2   Symmetry     Eigval 1'
1826     &      //'    Eigval 2     T1-norm        T2-norm'
1827
1828         ELSE
1829           WRITE (LUPRI,*) 'Unknown list in CC_TSTLST... ignored.'
1830           RETURN
1831         ENDIF
1832C
1833         DO ILSTNR = 1, NVEC
1834
1835           ISYM = ILSTSYM(LIST,ILSTNR)
1836           KT1  = 1
1837           KEND1= KT1  + NT1AM(ISYM)
1838Chol       IF (.NOT. CCS) THEN
1839           IF (.NOT. CCSEFF) THEN
1840              KT2   = KEND1
1841              KEND1 = KT2  + NT2AM(ISYM)
1842           ENDIF
1843           IF (CCR12) THEN
1844              KT12  = KEND1
1845              KEND1 = KT12  + NTR12AM(ISYM)
1846           END IF
1847           LEND1= LWORK - KEND1
1848           IF (LEND1 .LE. 0) THEN
1849               WRITE (LUPRI,*) 'Insufficient workspace in CC_TSTLST '
1850               CALL QUIT('Insufficient workspace in CC_TSTLST ')
1851           END IF
1852
1853           IOPT  = 3
1854           CALL CC_RDRSP(LIST,ILSTNR,ISYM,IOPT,MODEL,WORK(KT1),
1855     &                   WORK(KT2))
1856
1857           X1 = DDOT(NT1AM(ISYM),WORK(KT1),1,WORK(KT1),1)
1858Chol       IF (CCS) THEN
1859           IF (CCSEFF) THEN
1860              X2 = 0.0D0
1861           ELSE
1862              X2 = DDOT(NT2AM(ISYM),WORK(KT2),1,WORK(KT2),1)
1863           ENDIF
1864
1865           IF (CCR12) THEN
1866            IOPT = 32
1867            CALL CC_RDRSP(LIST,ILSTNR,ISYM,IOPT,MODEL,DUMMY,WORK(KT12))
1868            X12 = DDOT(NTR12AM(ISYM),WORK(KT12),1,WORK(KT12),1)
1869           ELSE
1870            X12 = 0.0D0
1871           END IF
1872
1873           IF     (LIST(2:2).NE.'C'.AND.ORDER.EQ.1.AND.NSTAT.EQ.0) THEN
1874
1875             ! ground state freq.-dep. first order response
1876             WRITE (LUPRI,'(1X,I5,5X,A8,2X,I4,4X,4(2X,F15.10))')
1877     &            ILSTNR,LABEL(ILSTNR,1),ISYM,FREQS(ILSTNR,1),X1,X2,X12
1878
1879           ELSEIF (LIST(1:1).NE.'C'.AND.ORDER.EQ.2.AND.NSTAT.EQ.0) THEN
1880
1881             ! ground state freq.-dep. second order response
1882             WRITE (LUPRI,'(1X,I5,5X,2(A8,2X),I4,6X,2F12.7,2X,
1883     &             F15.10,F15.10)')
1884     &            ILSTNR,LABEL(ILSTNR,1),LABEL(ILSTNR,2),ISYM,
1885     &                    FREQS(ILSTNR,1),FREQS(ILSTNR,2),X1,X2
1886
1887           ELSEIF (LIST(2:2).EQ.'C'.AND.ORDER.EQ.1.AND.NSTAT.EQ.0) THEN
1888
1889             ! cauchy expansion of ground state first order response
1890             WRITE (LUPRI,'(1X,I5,5X,A8,2X,I4,6X,I6,2X,E20.10,E20.10)')
1891     &            ILSTNR,LABEL(ILSTNR,1),ISYM,ICAU(ILSTNR,1),X1,X2
1892
1893           ELSEIF (LIST(1:1).EQ.'C'.AND.ORDER.EQ.2.AND.NSTAT.EQ.0) THEN
1894
1895             ! cauchy expansion of ground state second order response
1896             WRITE (LUPRI,'(1X,I5,5X,2(A8,2X),I4,2I10,2X,E20.10,
1897     &             E20.10)')
1898     &            ILSTNR,LABEL(ILSTNR,1),LABEL(ILSTNR,2),ISYM,
1899     &                     ICAU(ILSTNR,1), ICAU(ILSTNR,2),X1,X2
1900
1901           ELSEIF (ORDER.EQ.0 .AND. NSTAT.EQ.1) THEN
1902
1903            ! excited state vectors, no response
1904            IF (LIST(1:3).EQ.'E0 '.OR.LIST(1:3).EQ.'BE ') THEN
1905             WRITE (LUPRI,'(1X,I5,5X,I4,4X,I4,6X,F15.10,F15.10,F15.10)')
1906     &            ILSTNR,ISTAT(ILSTNR,1),ISYM,
1907     &            EIGVAL(ISTAT(ILSTNR,1),1),X1,X2
1908            ELSE
1909             WRITE (LUPRI,'(1X,I5,5X,I4,4X,I4,6X,F15.10,F15.10,F15.10)')
1910     &            ILSTNR,ISTAT(ILSTNR,1),ISYM,EIGVAL(ILSTNR,1),X1,X2
1911            END IF
1912
1913           ELSEIF (ORDER.EQ.1 .AND. NSTAT.EQ.1) THEN
1914
1915             ! excited state first order response
1916             WRITE (LUPRI,'(1X,I5,5X,A8,2X,I4,4X,I4,6X,4F15.10)')
1917     &            ILSTNR,LABEL(ILSTNR,1),ISTAT(ILSTNR,1),ISYM,
1918     &                  EIGVAL(ILSTNR,1),FREQS(ILSTNR,1),X1,X2
1919
1920           ELSEIF (ORDER.EQ.0 .AND. NSTAT.EQ.2) THEN
1921
1922             ! excited state - excited state vectors, no response
1923             WRITE (LUPRI,'(1X,I5,5X,2(I4,4X),I4,6X,2F15.10,F15.10,
1924     &             F15.10)')
1925     &            ILSTNR,ISTAT(ILSTNR,1),ISTAT(ILSTNR,2),ISYM,
1926     &                   EIGVAL(ILSTNR,1),EIGVAL(ILSTNR,2),X1,X2
1927
1928           ENDIF
1929
1930         END DO
1931
1932         IF (CHOINT .AND. CC2) THEN
1933            WRITE(LUPRI,'(/,1X,A,A)')
1934     &      'NOTICE: Cholesky CC2 doubles not stored! Artificial',
1935     &      ' zero norm printed.'
1936         ENDIF
1937
1938         WRITE (LUPRI,'(//)')
1939
1940         RETURN
1941         END
1942*=====================================================================*
1943