1
2* $Id$
3*
4
5
6* $Log: not supported by cvs2svn $
7* Revision 1.98  2009/03/19 20:42:17  bylaska
8* ...EJB
9*
10* Revision 1.97  2009/02/24 21:30:17  bert
11* In psi_finalize, occ1 and occ2 were deallocated in wrong order.
12*
13* Revision 1.96  2009/02/07 03:50:56  bylaska
14* Bassi Vectorization Fix...EJB
15*
16* Revision 1.95  2008/12/18 21:15:51  bylaska
17* ...updates for calculating spin contaminatio....EJB
18*
19* Revision 1.94  2008/11/17 17:25:45  bylaska
20* fractional occupation updates....EJB
21*
22* Revision 1.93  2008/10/22 23:56:43  bylaska
23* added NWCHEM_NWPW_LIBRARY to nwchemrc. fixed bug in paw...EJB
24*
25* Revision 1.92  2008/09/30 19:53:35  bylaska
26* Added Baden's exchange algorithm...EJB
27*
28* Revision 1.91  2008/09/17 00:55:36  bylaska
29* ...EJB
30*
31* Revision 1.90  2008/09/15 20:25:33  bylaska
32* ...fractional bug fix..EJB
33*
34* Revision 1.89  2008/09/11 21:26:51  bylaska
35* ...EJB
36*
37* Revision 1.88  2008/06/21 19:37:16  bylaska
38*
39* initalization error fixed with psi_Tgradient...EJB
40*
41* Revision 1.87  2008/06/02 15:20:04  bylaska
42* ..io fixes...EJB
43*
44* Revision 1.86  2008/05/13 02:10:36  bylaska
45* ...EJB
46*
47* Revision 1.85  2008/04/21 19:34:27  bylaska
48* queue fft added to cpsi_H, bug fixes in DMatrix_dgemm1_rot (MPI_Allgather routine replaced with a MPI_AllReduce for stability), np_orbital keyword replaced with np_dimensions keyword (needed for Parallel3d routines)
49* ...EJB
50*
51* Revision 1.84  2007/11/17 00:25:26  d3p708
52* ....PNJ
53*
54* Revision 1.83  2007/11/17 00:19:09  d3p708
55* pjn
56*    bugfix
57*
58* Revision 1.82  2007/11/16 22:32:53  d3p708
59* pjn
60*    stuff for berry phase dipole
61*
62* Revision 1.81  2007/10/01 23:02:46  bylaska
63* PAW changes..EJB
64*
65* Revision 1.80  2007/09/29 00:33:37  bylaska
66* ...EJB
67*
68* Revision 1.79  2007/09/24 16:58:14  bylaska
69* ...preliminary PAW modifications...
70*   - basis file format changed
71*   - .vpp formatting routines added to pspw
72*
73* - zdotc routines currently modified to tzdotc.
74* ...EJB
75*
76* Revision 1.78  2007/09/13 20:38:36  bylaska
77* occupation template added to band...EJB
78*
79* Revision 1.77  2007/03/27 02:02:49  bylaska
80* more qmmm_updates....EJB
81*
82* Revision 1.76  2007/03/22 20:46:22  bylaska
83* New implementation of QM/MM.
84* ....EJB
85*
86* Revision 1.75  2007/02/23 01:24:32  bylaska
87* ...EJB
88*
89* Revision 1.74  2007/02/10 03:56:54  bylaska
90* ...bug fix...
91* ..EJB
92*
93* Revision 1.73  2007/02/10 03:40:18  bylaska
94* replaced calls to Grsm_g_MakeOrtho with Dneall_f_ortho
95* ...EJB
96*
97* Revision 1.72  2007/01/02 18:36:52  bylaska
98* HGH pseudopotentials added to band.
99* ...EJB
100*
101* Revision 1.71  2006/10/13 01:43:58  bylaska
102* tcgmsg code for 2d grid distribution.  Also cleaned up Dmatrix_ calls so
103* that they are Dneall_ calls instead.
104* ....EJB
105*
106* Revision 1.70  2006/10/07 00:10:07  bylaska
107* Initial implementation of 2d processor grid parallelization in pspw.  Currently works with:
108*
109* task pspw steepest_descent
110* task pspw energy            (only minimizer 1, minimizer 2?, other minimizers not yet implemented)
111*
112* Currently only works with USE_MPIF option, straight tcgmsg only partially implemented.  Car-Parrinello, HFX, SIC, and various analysis codes are also not yet ported.
113*
114*
115* The number of processors along the orbital dimension, np_orbital, is entered as follows, e.g.:
116*
117* nwpw
118*    np_orbital  2
119* end
120*
121* The number of processors along the grid dimension, np_grid, is currently defined using np_orbital as
122*
123* np_grid = np/np_orbital
124*
125* where np is the total number of processors.
126*
127* ...EJB
128*
129* Revision 1.69  2006/09/20 19:18:49  bylaska
130* Adding Dmatrix
131* ...EJB
132*
133* Revision 1.68  2006/08/13 01:03:28  bylaska
134* Checking in code not include in 5.0 release.
135* A chain algorithm was added to Nose-Hoover thermostats.
136* Preliminary implementation of a processor group decomposition added to pspw, i.e. parallel decomposition is over fft grid and electrons.
137* ...EJB
138*
139* Revision 1.67  2006/01/26 18:29:36  bylaska
140* bug fix for gga checking...EJB
141*
142* Revision 1.66  2006/01/06 22:52:28  bylaska
143* parallel io bug fix...EJB
144*
145* Revision 1.65  2006/01/06 21:48:43  bylaska
146* io changes for inversion symmetry....EJB
147*
148* Revision 1.64  2005/12/29 03:06:09  marat
149* qmmm interface stuff
150*
151* Revision 1.63  2005/12/22 01:35:07  bylaska
152* revPBE added and gga logic restructured....EJB
153*
154* Revision 1.62  2005/10/05 21:21:30  bylaska
155* psi_iptr_write added...EJB
156*
157* Revision 1.61  2005/05/24 17:36:27  bylaska
158* Stresses added to SIC and HFX.
159* ....BLYP functional updates
160* ....ECCE hacks
161* ...EJB
162*
163* Revision 1.60  2005/02/09 02:39:10  bylaska
164* .............................EJB
165*
166* Revision 1.59  2005/01/17 20:51:33  edo
167* fixed a  couple of FPEs
168*
169* Revision 1.58.2.1  2005/01/17 20:51:06  edo
170* fixed a couple of FPEs
171*
172* Revision 1.58  2004/12/21 16:58:35  bylaska
173* various io fixes for dos...EJB
174*
175* Revision 1.57  2004/12/06 20:03:25  bylaska
176* RMM-DIIS diagonalizer added to pspw.
177* nwpw_list updated to handle multiple lists.
178* ....EJB
179*
180* Revision 1.56  2004/11/29 16:05:21  bylaska
181* Finite difference stresses added to PSPW, BAND, and PAW modules.
182*    - This is currently the default for BAND and PAW
183* Fixed the analytic unrestricted gga stress term in PSPW.
184* Fixed unrestricted optimization for minimizers 1 and 2 in the PSPW and PAW modules.
185* Partial implementation of analytic stress in BAND module.
186*    - kinetic, ewald, and coulomb stresses have been implemented
187*
188* ....EJB
189*
190* Revision 1.55  2004/11/08 01:32:49  bylaska
191*
192* Unrestricted and closed-shell restricted Hartree-Fock has been implemented into pspw.
193*    - works with minimizers 1,2,3,4,6,and 7.
194*    - band-by-band minimizer (minimizer 5) not yet implemented
195*    - free-space coulomb and screened-coulomb kernels implemented.
196*    - the free-space coulomb kernel has been tested.
197*    - the screened-coulomb kernel needs to be debugged/tested?
198*    - restricted open-shell HF has not yet been implemented
199*
200* ....EJB
201*
202* Revision 1.54  2004/03/08 22:51:27  bylaska
203* Fractional occupation working in pspw with minimizer 4, steepest descent, and Car-Parrinello.
204*
205* Bug fix in velocity initialization in liquid and solid-state Car-Parrinello simulations...incell3 instead of incell2 was used in newton step.
206*
207* Added restart capabilites to thermostat masses...Qe and Qr and eke0 are now propagated to
208* restart Car-Parrinello simulations.
209*
210* SIC input modifications.
211*
212* Wannier orbital output modifications.
213*
214* ....EJB
215*
216* Revision 1.53  2004/03/02 00:10:22  bylaska
217* ....EJB
218*
219* Revision 1.52  2004/03/01 05:14:33  bylaska
220* Mulliken and DOS fixes.
221* Added Mulliken projections based on atomic orbitals
222* Added projected density of states (based on Mulliken projections)
223* ...EJB
224*
225* Revision 1.51  2004/02/06 01:57:22  bylaska
226*
227* Option added to write out temporary psi for Kiril.
228* Tempory psi written if
229* set nwpw:psi_tmp .true.
230*
231*
232* ....EJB
233*
234* Revision 1.50  2004/01/28 00:08:47  bylaska
235* Bug fixes...EJB
236*
237* Revision 1.49  2003/12/13 21:07:41  bylaska
238* Kohn-Sham scf minimizer added to pspw. Mixing includes
239* 	- simple mixing
240* 	- Anderson mixing
241* 	- Johnson-Pulay mixing
242* 	- Kerker preconditioner
243*
244* ....EJB
245*
246* Revision 1.48  2003/12/02 23:55:34  bylaska
247* density of state generation added...EJB
248*
249* Revision 1.47  2003/12/02 19:17:10  bylaska
250* HGH pseudpotential added.
251* TM, Hamman, HGH, pspw_default, and paw_default pseudopotential libraries have been added.
252* KS minimizer updates.
253* ...EJB
254*
255
256
257************************ f_orb orbitals Part ************************
258
259*     ***********************************
260*     *				        *
261*     *	     psi_minimize_f_orb         *
262*     *				        *
263*     ***********************************
264      subroutine psi_minimize_f_orb()
265      implicit none
266
267#include "bafdecls.fh"
268#include "psi.fh"
269
270      !*** local variables ***
271      integer maxit_orb
272      integer ii,l,l2
273      real*8  sum,maxerror,error_out,e0
274
275      !*** external functions ***
276      real*8   control_tole
277      external control_tole
278
279      !call psi_gen_density_potentials(1)
280      maxit_orb=120
281      maxerror = control_tole()
282
283      do ii=1,(ne(1)+ne(2))
284         l2 = 0
285
286         !*** orthogonalize to lower orbitals  ****
287 2       l2 = l2 + 1
288         call psi_project_out_f_orb1(
289     >           ii,
290     >           dcpl_mb(psi1(1)+(ii-1)*npack1))
291
292         !*** normalize ****
293         call Pack_cc_dot(1,
294     >            dcpl_mb(psi1(1) +(ii-1)*npack1),
295     >            dcpl_mb(psi1(1) +(ii-1)*npack1),
296     >            sum)
297         sum = 1.0d0/dsqrt(sum)
298         call Pack_c_SMul1(1,sum,dcpl_mb(psi1(1) +(ii-1)*npack1))
299
300         !*** minimize orbital ****
301          l = 0
302 3        call psi_KS_update_f_orb(maxit_orb,
303     >                               maxerror,
304     >                               0.001d0,ii,error_out,e0)
305          !write(*,*) "e0:",ii,l,e0,error_out
306          l = l+1
307          if ((error_out.gt.maxerror).and.(l.le.(1+(l2-1)*3))) go to 3
308          if (((error_out.gt.maxerror).or.(e0.gt.4.0d0))
309     >        .and.(l2.le.1)) then
310           call Pack_c_Zero(1,dcpl_mb(psi1(1) +(ii-1)*npack1))
311           call Pack_c_setzero(1,1.0d0,dcpl_mb(psi1(1) +(ii-1)*npack1))
312           go to 2
313          end if
314
315          dbl_mb(eig(1)+ii-1) = e0
316
317      end do
318      call psi_sort_f_orb()
319
320
321      return
322      end
323
324      subroutine psi_sort_f_orb()
325      implicit none
326#include "errquit.fh"
327
328#include "bafdecls.fh"
329#include "psi.fh"
330
331      logical value
332      integer i,j,ii,jj,ms
333      integer r1(2)
334      real*8  ei,ej
335
336      value = BA_push_get(mt_dcpl,npack1,'r1',r1(2),r1(1))
337      if (.not. value) call errquit(
338     >     'psi_sort_f_orb: out of stack memory',0,MA_ERR)
339
340      do ms=1,ispin
341
342        !*** Bubble sort ***
343        do ii=1,ne(ms)
344         do jj=ii+1,ne(ms)
345           i = ii + (ms-1)*ne(1)
346           j = jj + (ms-1)*ne(1)
347           ei = dbl_mb(eig(1)+i-1)
348           ej = dbl_mb(eig(1)+j-1)
349
350           !*** swap ***
351           if (ej.lt.ei) then
352             dbl_mb(eig(1)+i-1) = ej
353             dbl_mb(eig(1)+j-1) = ei
354             call Pack_c_Copy(1,dcpl_mb(psi1(1)+(i-1)*npack1),
355     >                          dcpl_mb(r1(1)))
356             call Pack_c_Copy(1,dcpl_mb(psi1(1)+(j-1)*npack1),
357     >                          dcpl_mb(psi1(1)+(i-1)*npack1))
358             call Pack_c_Copy(1,dcpl_mb(r1(1)),
359     >                          dcpl_mb(psi1(1)+(j-1)*npack1))
360           end if
361
362         end do
363        end do
364
365      end do
366
367      value = BA_pop_stack(r1(2))
368      if (.not. value) call errquit(
369     >     'psi_sort_f_orb: popping stack memory',1, MA_ERR)
370      return
371      end
372
373*     ***********************************
374*     *				        *
375*     *	     psi_KS_update_f_orb      *
376*     *				        *
377*     ***********************************
378
379*    This routine performs a KS update on orbital ii
380*
381      subroutine psi_KS_update_f_orb(maxiteration,
382     >                             maxerror,perror,ii,
383     >                             error_out,e0)
384      implicit none
385#include "errquit.fh"
386      integer maxiteration
387      real*8  maxerror,perror
388      integer ii
389      real*8 error_out
390      real*8 e0
391
392#include "bafdecls.fh"
393#include "psi.fh"
394
395*     **** local variables ****
396      integer MASTER,taskid
397      parameter (MASTER=0)
398
399      logical value,done,oneloop
400      integer it
401      real*8 eold,percent_error,error0,de0,lmbda_r0,lmbda_r1
402      real*8 theta
403      integer r1(2),t0(2),t(2),g(2)
404      integer psi_ptr
405
406      psi_ptr=psi1(1)
407
408      call Parallel_taskid(taskid)
409
410      value = BA_push_get(mt_dcpl,npack1,'t0',t0(2),t0(1))
411      value = value.and.
412     >        BA_push_get(mt_dcpl,npack1,'r1',r1(2),r1(1))
413      value = value.and.
414     >        BA_push_get(mt_dcpl,npack1,'g',g(2),g(1))
415      value = value.and.
416     >        BA_push_get(mt_dcpl,npack1,'t',t(2),t(1))
417      if (.not. value) call errquit(
418     >     'psi_KS_update_f_orb: out of stack memory',0, MA_ERR)
419
420      done = .false.
421      error0 = 0.0d0
422      e0 = 0.0d0
423      theta = -3.14159d0/600.0d0
424      lmbda_r0 = 1.0d0
425      it = 0
426 2    continue
427
428         it = it + 1
429         eold = e0
430
431*        *** calculate residual (steepest descent) direction for a single band ***
432         call psi_get_gradient_f_orb(ii,dcpl_mb(g(1)))
433         call Pack_cc_dot(1,dcpl_mb(psi_ptr+(ii-1)*npack1),
434     >                   dcpl_mb(g(1)),
435     >                    e0)
436
437         e0 = -e0
438         percent_error=0d0
439         if(error0.ne.0d0)
440     A      percent_error = dabs(e0-eold)/error0
441
442         done = ((it.gt.maxiteration)
443     >           .or.
444     >           (dabs(e0-eold).lt.maxerror))
445
446         if (done) go to 4
447
448
449         call Pack_c_Copy(1,dcpl_mb(g(1)),dcpl_mb(r1(1)))
450         call Pack_cc_daxpy(1,(e0),
451     >                 dcpl_mb(psi_ptr+(ii-1)*npack1),
452     >                 dcpl_mb(r1(1)))
453
454
455         !*** determine conjuagate direction ***
456         call Pack_cc_dot(1,dcpl_mb(r1(1)),
457     >                   dcpl_mb(r1(1)),
458     >                   lmbda_r1)
459         call Pack_c_Copy(1,dcpl_mb(r1(1)),dcpl_mb(t(1)))
460
461         if (it.gt.1) then
462         call Pack_cc_daxpy(1,(lmbda_r1/lmbda_r0),
463     >                   dcpl_mb(t0(1)),
464     >                   dcpl_mb(t(1)))
465         end if
466         lmbda_r0 = lmbda_r1
467         oneloop = .true.
468 3       call Pack_c_Copy(1,dcpl_mb(t(1)),dcpl_mb(t0(1)))
469
470
471*        *** normalize search direction, t ****
472         call psi_project_out_f_orb(ii,dcpl_mb(t(1)))
473         call Pack_cc_dot(1,dcpl_mb(t(1)),
474     >                   dcpl_mb(t(1)),
475     >                   de0)
476         de0 = 1.0d0/dsqrt(de0)
477c         call Pack_c_SMul(1,de0,dcpl_mb(t(1)),dcpl_mb(t(1)))
478         call Pack_c_SMul1(1,de0,dcpl_mb(t(1)))
479         call Pack_cc_dot(1,dcpl_mb(t(1)),
480     >                   dcpl_mb(g(1)),
481     >                   de0)
482
483*        *** bad direction ***
484         if ((de0.lt.0.0d0).and.oneloop) then
485           call Pack_c_Copy(1,dcpl_mb(g(1)),dcpl_mb(t(1)))
486           oneloop = .false.
487           go to 3
488         end if
489
490         de0 = -2.0d0*de0
491         call psi_linesearch_f_orb(ii,
492     >                               theta,e0,de0,dcpl_mb(t(1)))
493
494      go to 2
495
496
497*     **** release stack memory ****
498 4    value =           BA_pop_stack(t(2))
499      value = value.and.BA_pop_stack(g(2))
500      value = value.and.BA_pop_stack(r1(2))
501      value = value.and.BA_pop_stack(t0(2))
502      if (.not. value) call errquit(
503     >     'psi_KS_update_f_orb: popping stack memory',1, MA_ERR)
504
505      error_out = dabs(e0-eold)
506      e0        = -e0
507      return
508      end
509
510*     ***********************************
511*     *				        *
512*     *	     psi_linesearch_f_orb
513*     *				        *
514*     ***********************************
515
516*    This routine performs a linesearch on orbital ii, in the direction t.
517* This routine is needed for a KS minimizer.
518*  e0 = <orb|g>
519*  de0 = 2*<t|g>
520*
521      subroutine psi_linesearch_f_orb(ii,theta,e0,de0,t)
522      implicit none
523#include "errquit.fh"
524      integer ii
525      real*8  theta
526      real*8  e0,de0
527      complex*16 t(*) !search direction
528
529#include "bafdecls.fh"
530#include "psi.fh"
531
532*     **** local variables ****
533      logical value
534      integer orb(2),g(2),psi_ptr
535      real*8 x,y,pi,dtheta_min,e1
536
537      psi_ptr=psi1(1)
538
539      pi = 4.0d0*datan(1.0d0)
540      !dtheta = pi/300.0d0
541      dtheta_min = 0.01*theta
542
543*     **** allocate stack memory ****
544      value = BA_push_get(mt_dcpl,npack1,'orb',
545     >                       orb(2),orb(1))
546      value = value.and.
547     >        BA_push_get(mt_dcpl,npack1,'g',
548     >                       g(2),g(1))
549      if (.not. value) call errquit(
550     >     'psi_linesearch_f_orb: out of stack memory',0, MA_ERR)
551
552
553      call Pack_c_Copy(1,dcpl_mb(psi_ptr+(ii-1)*npack1),
554     >                   dcpl_mb(orb(1)))
555
556*     **** orb2 = orb*cos(pi/300) + t*sin(pi/300) ****
557  10  x = cos(theta)
558      y = sin(theta)
559      call Pack_c_SMul(1,x,
560     >                  dcpl_mb(orb(1)),
561     >                  dcpl_mb(psi_ptr+(ii-1)*npack1))
562      call Pack_cc_daxpy(1,y,
563     >                   t,
564     >                   dcpl_mb(psi_ptr+(ii-1)*npack1))
565
566*     *** determine theta ***
567      call psi_get_gradient_f_orb(ii,dcpl_mb(g(1)))
568      call Pack_cc_dot(1,dcpl_mb(psi_ptr+(ii-1)*npack1),
569     >                   dcpl_mb(g(1)),
570     >                   e1)
571      e1 = -e1
572
573
574
575      x = (e0 - e1 + 0.5d0*de0*sin(2*theta))
576     >    /(1.0d0-cos(2*theta))
577      theta = 0.5d0*datan(0.5d0*de0/x)
578
579
580
581*     **** orb2 = orb*cos(theta) + t*sin(theta) ****
582      x = cos(theta)
583      y = sin(theta)
584      call Pack_c_SMul(1,x,
585     >                  dcpl_mb(orb(1)),
586     >                  dcpl_mb(psi_ptr+(ii-1)*npack1))
587      call Pack_cc_daxpy(1,y,
588     >                   t,
589     >                   dcpl_mb(psi_ptr+(ii-1)*npack1))
590
591
592*     **** release stack memory ****
593      value =           BA_pop_stack(g(2))
594      value = value.and.BA_pop_stack(orb(2))
595      if (.not. value) call errquit(
596     >     'psi_linesearch_f_orb: popping stack memory',1,MA_ERR)
597
598      return
599      end
600
601
602*     ***********************************
603*     *				        *
604*     *	     psi_get_gradient_f_orb	*
605*     *				        *
606*     ***********************************
607
608*    This routine returns the Hpsi(i).
609* This routine is needed for a KS minimizer.
610*
611      subroutine psi_get_gradient_f_orb(ii,Horb)
612      implicit none
613      integer ii
614      complex*16 Horb(*)
615
616#include "bafdecls.fh"
617#include "psi.fh"
618
619*     **** local variables ****
620      integer psi_ptr,ms
621
622      psi_ptr=psi1(1)+(ii-1)*npack1
623
624      if (ii.le.neq(1)) then
625         ms = 1
626      else
627         ms = 2
628      end if
629
630      call electron_get_gradient_virtual(ms,dcpl_mb(psi_ptr),Horb)
631
632      return
633      end
634
635*     *******************************************
636*     *				                *
637*     *	         psi_project_out_f_orb1        *
638*     *				                *
639*     *******************************************
640*
641*    This routine projects out non-orthogonal components of Horb.
642* This routine is needed for a KS minimizer.
643*
644      subroutine psi_project_out_f_orb1(ii,Horb)
645      integer ii
646      complex*16 Horb(*)
647
648#include "bafdecls.fh"
649#include "psi.fh"
650
651      integer ms,jj,kk,shift,shifte
652      real*8  sum
653
654*     **** spin up orbital ****
655      if (ii.le.neq(1)) then
656
657         shift  = 0
658         shifte = 0
659         ms     = 1
660         kk     = ii
661*     **** spin down orbital ****
662      else
663         shift  = neq(1)*npack1
664         shifte = neq(1)*npack1
665         ms     = 2
666         kk     = ii-neq(1)
667      end if
668
669
670      !**** project out orbitals ****
671      do jj=1,(kk-1)
672        call Pack_cc_dot(1,
673     >            dcpl_mb(psi1(1) +(jj-1)*npack1+shifte),
674     >            Horb,
675     >            sum)
676
677        call Pack_cc_daxpy(1,(-sum),
678     >            dcpl_mb(psi1(1) +(jj-1)*npack1+shifte),
679     >            Horb)
680      end do
681
682
683      return
684      end
685
686
687
688
689*     *******************************************
690*     *				                *
691*     *	         psi_project_out_f_orb          *
692*     *				                *
693*     *******************************************
694*
695*    This routine projects out non-orthogonal components of Horb.
696* This routine is needed for a KS minimizer.
697*
698      subroutine psi_project_out_f_orb(ii,Horb)
699      integer ii
700      complex*16 Horb(*)
701
702#include "bafdecls.fh"
703#include "psi.fh"
704
705      integer ms,jj,kk,shift,shifte
706      real*8  sum
707
708*     **** spin up orbital ****
709      if (ii.le.neq(1)) then
710
711         shift  = 0
712         shifte = 0
713         ms     = 1
714         kk     = ii
715*     **** spin down orbital ****
716      else
717         shift  = neq(1)*npack1
718         shifte = neq(1)*npack1
719         ms     = 2
720         kk     = ii-neq(1)
721      end if
722
723
724      !**** project out  orbitals ****
725      do jj=1,(kk)
726        call Pack_cc_dot(1,
727     >            dcpl_mb(psi1(1) +(jj-1)*npack1+shifte),
728     >            Horb,
729     >            sum)
730
731        call Pack_cc_daxpy(1,(-sum),
732     >            dcpl_mb(psi1(1) +(jj-1)*npack1+shifte),
733     >            Horb)
734      end do
735
736
737      return
738      end
739
740
741
742************************ gradient virtural orbital Part ************************
743*     ***********************************
744*     *                                 *
745*     *      psi_gen_hml_virtual        *
746*     *                                 *
747*     ***********************************
748      subroutine psi_gen_hml_virtual(assending)
749      implicit none
750      logical assending
751
752#include "bafdecls.fh"
753#include "psi.fh"
754#include "errquit.fh"
755
756      !*** local variables ***
757      integer ii,jj,jjstart,jjend,ms,mshift,indx,indxt
758      integer Horb(2)
759      real*8  tsum
760
761*     **** allocate temporary memory from stack ****
762      if (.not.BA_push_get(mt_dcpl,npack1,'Horb',Horb(2),Horb(1)))
763     >   call errquit('psi_gen_hml_virtual:pusing stack',0,MA_ERR)
764
765      do ii=1,(ne_excited(1)+ne_excited(2))
766
767         !*** orthogonalize to lower orbitals  ****
768         call psi_project_out_virtual1(
769     >           ii,
770     >           dcpl_mb(psi1_excited(1)+(ii-1)*npack1))
771
772         !*** normalize ****
773         call Pack_cc_dot(1,
774     >            dcpl_mb(psi1_excited(1) +(ii-1)*npack1),
775     >            dcpl_mb(psi1_excited(1) +(ii-1)*npack1),
776     >            tsum)
777         tsum = 1.0d0/dsqrt(tsum)
778         call Pack_c_SMul1(1,tsum,
779     >            dcpl_mb(psi1_excited(1) +(ii-1)*npack1))
780
781         if (ii.le.ne_excited(1)) then
782            ms = 1
783            jjstart = 1
784            jjend   = ii
785            mshift = 0
786         else
787            ms = 2
788            jjstart = ne_excited(1)+1
789            jjend   = ii
790            mshift = ne_excited(1)*ne_excited(1)
791         end if
792c         call Pack_cc_dot(1,
793c     >               dcpl_mb(psi1_excited(1) +(ii-1)*npack1),
794c     >               dcpl_mb(psi1_excited(1) +(ii-1)*npack1),
795c     >               tsum)
796         !write(*,*) "ii,norm=",ii,tsum
797
798         call electron_get_gradient_virtual(ms,
799     >             dcpl_mb(psi1_excited(1)+(ii-1)*npack1),
800     >             dcpl_mb(Horb(1)))
801
802         do jj=jjstart,jjend
803c         call Pack_cc_dot(1,
804c     >               dcpl_mb(psi1_excited(1) +(ii-1)*npack1),
805c     >               dcpl_mb(psi1_excited(1) +(jj-1)*npack1),
806c     >               tsum)
807         !write(*,*) "ii,jj,norm=",ii,jj,tsum
808            call Pack_cc_dot(1,
809     >               dcpl_mb(psi1_excited(1) +(jj-1)*npack1),
810     >               dcpl_mb(Horb(1)),
811     >               tsum)
812            !write(*,*) "ii,jj,e=",ii,jj,tsum
813            indx  = (ii-jjstart)+(jj-jjstart)*ne_excited(ms)+mshift
814            dbl_mb(hml_excited(1)+indx)  = tsum
815            if (ii.ne.jj) then
816               indxt = (jj-jjstart)+(ii-jjstart)*ne_excited(ms)+mshift
817               dbl_mb(hml_excited(1)+indxt) = tsum
818            end if
819         end do
820      end do
821
822*     **** deallocate temporary memory from stack ****
823      if (.not.BA_pop_stack(Horb(2)))
824     >   call errquit('psi_gen_hml_virtual:popping stack',0,MA_ERR)
825
826c      call Dnexall_m_diagonalize(0,ispin,ne_excited,
827c     >                           dbl_mb(hml_excited(1)),
828c     >                           dbl_mb(eig_excited(1)),.true.)
829      call Dnexall_m_diagonalize(0,ispin,ne_excited,
830     >                           dbl_mb(hml_excited(1)),
831     >                           dbl_mb(eig_excited(1)),assending)
832      call  Dnexall_fmf_Multiply(0,ispin,ne_excited,
833     >                           dcpl_mb(psi1_excited(1)),npack1,
834     >                           dbl_mb(hml_excited(1)),1.0d0,
835     >                           dcpl_mb(psi2_excited(1)),0.0)
836      indx = psi1_excited(1)
837      psi1_excited(1) = psi2_excited(1)
838      psi2_excited(1) = indx
839c      call dcopy(2*npack1,
840c     >           dcpl_mb(psi1_excited(1)+3*npack1),1,
841c     >           dcpl_mb(psi1(1)),1)
842
843
844      return
845      end
846
847
848
849************************ CI virtural orbital Part ************************
850*     ***********************************
851*     *                                 *
852*     *    psi_minimize_virtual_CI      *
853*     *                                 *
854*     ***********************************
855      subroutine psi_minimize_virtual_CI(ci_algorithm)
856      implicit none
857      integer ci_algorithm
858
859#include "bafdecls.fh"
860#include "errquit.fh"
861#include "util.fh"
862#include "stdio.fh"
863#include "psi.fh"
864
865      real*8 Hgg
866      common / CI_Hgg_common / Hgg
867
868      !*** local variables ***
869      integer taskid,MASTER
870      parameter (MASTER=0)
871
872      logical oprint
873      integer maxit_sweeps,maxit_orb,maxit_ls
874      integer ii,l
875      real*8  sum,maxerror,error_out,e0,eci
876      integer g93(2),eig_eci(2)
877
878      !*** external functions ***
879      logical  control_print
880      external control_print
881      integer  control_CI_maxit_orb,control_CI_maxit_linesearch
882      external control_CI_maxit_orb,control_CI_maxit_linesearch
883      integer  control_CI_maxit_sweeps
884      external control_CI_maxit_sweeps
885      real*8   control_tole
886      external control_tole
887
888      call Parallel_taskid(taskid)
889      oprint= ((taskid.eq.MASTER).and.control_print(print_medium))
890
891      Hgg = 0.0d0
892
893      if (.not.BA_push_get(mt_dcpl,npack1,'g93',g93(2),g93(1)))
894     >   call errquit('psi_CI_minimize_virtual:out of stack',0,MA_ERR)
895      if (.not.BA_push_get(mt_dbl,ne_excited(1)+ne_excited(2),
896     >                     'eig_eci',eig_eci(2),eig_eci(1)))
897     >   call errquit('psi_CI_minimize_virtual:out of stack',1,MA_ERR)
898
899      !call psi_gen_density_potentials(1)
900      maxit_sweeps = control_CI_maxit_sweeps()
901      maxit_orb = control_CI_maxit_orb()
902      maxit_ls  = control_CI_maxit_linesearch()
903      maxerror  = control_tole()
904      if (oprint) then
905         write(luout,*)
906         write(luout,*) "COVOs Mimimization"
907         write(luout,*) "------------------"
908         write(luout,'(A,I9)') "CI gradient algorithm = ",ci_algorithm
909         write(luout,'(A,I9)') "maxit_sweeps          = ",maxit_sweeps
910         write(luout,'(A,I9)') "maxit_orb             = ",maxit_orb
911         write(luout,'(A,I9)') "maxit_linesearch      = ",maxit_ls
912         write(luout,'(A,E9.3)') "maxerror              = ",maxerror
913         write(luout,*)
914      end if
915      !maxerror = 1.0d-12
916      eci = 9.99d9
917      e0  = 9.99d9
918
919      do ii=1,(ne_excited(1)+ne_excited(2))
920
921         !*** orthogonalize to lower orbitals  ****
922         call psi_project_out_virtual1(
923     >           ii,
924     >           dcpl_mb(psi1_excited(1)+(ii-1)*npack1))
925
926         !*** normalize ****
927         call Pack_cc_dot(1,
928     >            dcpl_mb(psi1_excited(1) +(ii-1)*npack1),
929     >            dcpl_mb(psi1_excited(1) +(ii-1)*npack1),
930     >            sum)
931         sum = 1.0d0/dsqrt(sum)
932         call Pack_c_SMul1(1,sum,
933     >            dcpl_mb(psi1_excited(1) +(ii-1)*npack1))
934
935         call psi_get_gradient_virtual(ii,dcpl_mb(g93(1)))
936         call Pack_cc_dot(1,dcpl_mb(psi1_excited(1)+(ii-1)*npack1),
937     >                      dcpl_mb(g93(1)),e0)
938         e0 = e0
939
940         !*** minimize orbital ****
941
942          l = 0
943 2        call psi_CI_update_virtual(ci_algorithm,
944     >                               maxit_orb,maxit_ls,
945     >                               maxerror,
946     >                               0.001d0,ii,error_out,e0,eci)
947          l = l+1
948          !if ((error_out.gt.maxerror).and.(l.le.4)) go to 2
949          if ((error_out.gt.maxerror).and.(l.lt.maxit_sweeps)) go to 2
950
951          dbl_mb(eig_excited(1)+ii-1) = e0
952          dbl_mb(eig_eci(1)+ii-1)     = eci
953
954      end do
955      call psi_sort_virtual_CI(dbl_mb(eig_eci(1)))
956
957c      if (oprint) then
958c         do ii=1,(ne_excited(1)+ne_excited(2))
959c            write(*,*) "ii,eig=",ii,dbl_mb(eig_eci(1)+ii-1),
960c     >                              dbl_mb(eig_excited(1)+ii-1)
961c         end do
962c      end if
963
964      if (.not.BA_pop_stack(eig_eci(2)))
965     >   call errquit('psi_CI_minimize_virtual:pop stack',1,MA_ERR)
966      if (.not.BA_pop_stack(g93(2)))
967     >   call errquit('psi_CI_minimize_virtual:pop stack',2,MA_ERR)
968
969
970      return
971      end
972
973      subroutine psi_sort_virtual_CI(eig_eci)
974      implicit none
975      real*8 eig_eci(*)
976
977#include "bafdecls.fh"
978#include "psi.fh"
979#include "errquit.fh"
980
981      logical value
982      integer i,j,ii,jj,ms
983      integer r1(2)
984      real*8  ei,ej,eci,ecj
985
986      value = BA_push_get(mt_dcpl,npack1,'rr1',r1(2),r1(1))
987      if (.not. value) call errquit(
988     >     'psi_sort_virtual: out of stack memory',0, MA_ERR)
989
990      do ms=1,ispin
991
992        !*** Bubble sort ***
993        do ii=1,ne_excited(ms)
994         do jj=ii+1,ne_excited(ms)
995           i = ii + (ms-1)*ne_excited(1)
996           j = jj + (ms-1)*ne_excited(1)
997           ei = dbl_mb(eig_excited(1)+i-1)
998           ej = dbl_mb(eig_excited(1)+j-1)
999           eci = eig_eci(i)
1000           ecj = eig_eci(j)
1001
1002           !*** swap ***
1003           if (ecj.lt.eci) then
1004             eig_eci(i) = ecj
1005             eig_eci(j) = eci
1006             dbl_mb(eig_excited(1)+i-1) = ej
1007             dbl_mb(eig_excited(1)+j-1) = ei
1008             call Pack_c_Copy(1,dcpl_mb(psi1_excited(1)+(i-1)*npack1),
1009     >                          dcpl_mb(r1(1)))
1010             call Pack_c_Copy(1,dcpl_mb(psi1_excited(1)+(j-1)*npack1),
1011     >                          dcpl_mb(psi1_excited(1)+(i-1)*npack1))
1012             call Pack_c_Copy(1,dcpl_mb(r1(1)),
1013     >                          dcpl_mb(psi1_excited(1)+(j-1)*npack1))
1014           end if
1015
1016         end do
1017        end do
1018
1019      end do
1020
1021      value = BA_pop_stack(r1(2))
1022      if (.not. value) call errquit(
1023     >     'psi_sort_virtual: popping stack memory',1, MA_ERR)
1024      return
1025      end
1026
1027
1028
1029*     ***********************************
1030*     *                                 *
1031*     *      psi_CI_update_virtual      *
1032*     *                                 *
1033*     ***********************************
1034
1035*    This routine performs a KS update on virtual ii
1036*
1037      subroutine psi_CI_update_virtual(ci_algorithm,maxiteration,
1038     >                             maxit_ls,
1039     >                             maxerror,perror,ii,
1040     >                             error_out,e0,eci)
1041      implicit none
1042      integer ci_algorithm
1043      integer maxiteration,maxit_ls
1044      real*8  maxerror,perror
1045      integer ii
1046      real*8 error_out
1047      real*8 e0,eci
1048
1049#include "bafdecls.fh"
1050#include "psi.fh"
1051#include "errquit.fh"
1052#include "util.fh"
1053#include "stdio.fh"
1054
1055*     **** local variables ****
1056      integer MASTER,taskid
1057      parameter (MASTER=0)
1058
1059      logical value,done,oneloop,precondition,oprint
1060      integer it,pit,n2ft3d,nx,ny,nz
1061      real*8 eold,percent_error,error0,de0,lmbda_r0,lmbda_r1
1062      real*8 theta,ep,sp,e0old
1063      integer r1(2),t0(2),t(2),g(2),psig_r(2),psie_r(2)
1064      integer rhog(2),rhoe(2),rhoge(2),vcgg(2),vcge(2),vcee(2)
1065      integer vc(2),vctmp(2),te(2)
1066      integer psig_ptr,psie_ptr,psig_homo(2),hpsig(2),hpsig_r(2)
1067      real*8  Sgg,Sge,See
1068      real*8  H1gg,H1ge,H1ee
1069      real*8  H2gg,H2ge,H2ee
1070      real*8  Hgg,Hge,Hee
1071      real*8  scal1,scal2,dv,ehartr
1072      real*8  A,B,C,Elow,Ehigh,c1,c2,lmbda,sum,sigma
1073
1074      logical  control_print
1075      external control_print
1076      real*8   control_Ep,control_Sp,lattice_omega
1077      external control_Ep,control_Sp,lattice_omega
1078      integer  electron_psi_r_ptr
1079      external electron_psi_r_ptr
1080
1081      psig_ptr=psi1(1)
1082      psie_ptr=psi1_excited(1) + (ii-1)*npack1
1083
1084      call Parallel_taskid(taskid)
1085      oprint= ((taskid.eq.MASTER).and.control_print(print_medium))
1086
1087      call D3dB_nx(1,nx)
1088      call D3dB_ny(1,ny)
1089      call D3dB_nz(1,nz)
1090      scal1 = 1.0d0/dble(nx*ny*nz)
1091      scal2 = 1.0d0/lattice_omega()
1092      dv = scal1*lattice_omega()
1093
1094      call D3dB_n2ft3d(1,n2ft3d)
1095
1096      value = BA_push_get(mt_dcpl,npack1,'t0',t0(2),t0(1))
1097      value = value.and.
1098     >        BA_push_get(mt_dcpl,npack1,'r1',r1(2),r1(1))
1099      value = value.and.
1100     >        BA_push_get(mt_dcpl,npack1,'g',g(2),g(1))
1101      value = value.and.
1102     >        BA_push_get(mt_dcpl,npack1,'t',t(2),t(1))
1103      value = value.and.
1104     >        BA_push_get(mt_dcpl,npack1,'te',te(2),te(1))
1105
1106      value = value.and.
1107     >        BA_push_get(mt_dbl,n2ft3d,'psig_r',psig_r(2),psig_r(1))
1108      value = value.and.
1109     >        BA_push_get(mt_dbl,n2ft3d,'psie_r',psie_r(2),psie_r(1))
1110      value = value.and.
1111     >        BA_push_get(mt_dbl,n2ft3d,'rhog',rhog(2),rhog(1))
1112      value = value.and.
1113     >        BA_push_get(mt_dbl,n2ft3d,'rhoe',rhoe(2),rhoe(1))
1114      value = value.and.
1115     >        BA_push_get(mt_dbl,n2ft3d,'rhoge',rhoge(2),rhoge(1))
1116
1117      value = value.and.
1118     >        BA_push_get(mt_dbl,n2ft3d,'vctmp',vctmp(2),vctmp(1))
1119      value = value.and.
1120     >        BA_push_get(mt_dbl,n2ft3d,'vc',vc(2),vc(1))
1121
1122      value = value.and.
1123     >        BA_push_get(mt_dbl,n2ft3d,'vcgg',vcgg(2),vcgg(1))
1124      value = value.and.
1125     >        BA_push_get(mt_dbl,n2ft3d,'vcge',vcge(2),vcge(1))
1126      value = value.and.
1127     >        BA_push_get(mt_dbl,n2ft3d,'vcee',vcee(2),vcee(1))
1128      if (ci_algorithm.eq.5) then
1129         value = value.and.
1130     >           BA_push_get(mt_dcpl,npack1,'psig_homo',
1131     >                       psig_homo(2),psig_homo(1))
1132         value = value.and.
1133     >           BA_push_get(mt_dcpl,(neq(1)+neq(2))*npack1,'hpsig',
1134     >                       hpsig(2),hpsig(1))
1135         value = value.and.
1136     >           BA_push_get(mt_dbl,(neq(1)+neq(2))*n2ft3d,'hpsig_r',
1137     >                       hpsig_r(2),hpsig_r(1))
1138         call Pack_c_Copy(1,dcpl_mb(psig_ptr),dcpl_mb(psig_homo(1)))
1139      end if
1140      if (.not.value) call errquit(
1141     >     'psi_CI_update_virtual: out of stack memory',0, MA_ERR)
1142
1143      call Pack_c_Copy(1,dcpl_mb(psig_ptr),dbl_mb(psig_r(1)))
1144      call Pack_c_unpack(1,dbl_mb(psig_r(1)))
1145      call D3dB_cr_fft3b(1,dbl_mb(psig_r(1)))
1146      call D3dB_r_Zero_Ends(1,dbl_mb(psig_r(1)))
1147
1148      ep = control_Ep()
1149      sp = control_Sp()
1150      !precondition = .true.
1151      precondition = .false.
1152      done = .false.
1153      error0 = 0.0d0
1154      eci = 0.0d0
1155      !theta = 3.14159d0/600.0d0
1156      theta = 3.14159d0/10.0d0
1157      lmbda_r0 = 1.0d0
1158      it = 0
1159      pit = 0
1160 2    continue
1161         it = it + 1
1162
1163         if (ci_algorithm.eq.2) then
1164*        *** calculate residual (steepest descent) direction for psie ***
1165         call psi_2x2_virtual_gradient(dcpl_mb(psig_ptr),
1166     >                                 dbl_mb(psig_r(1)),
1167     >                                 dcpl_mb(psie_ptr),
1168     >                                 dbl_mb(psie_r(1)),
1169     >                                 eci,dcpl_mb(g(1)),
1170     >                                 dcpl_mb(t(1)),
1171     >                                 dbl_mb(rhog(1)),
1172     >                                 dbl_mb(vcee(1)),dbl_mb(vcge(1)))
1173         else if (ci_algorithm.eq.3) then
1174         call psi_3x3_virtual_gradient(dcpl_mb(psig_ptr),
1175     >                                 dbl_mb(psig_r(1)),
1176     >                                 dcpl_mb(psie_ptr),
1177     >                                 dbl_mb(psie_r(1)),
1178     >                                 eci,dcpl_mb(g(1)),
1179     >                                 dcpl_mb(t(1)),
1180     >                                 dcpl_mb(te(1)),
1181     >                                 dbl_mb(rhog(1)),
1182     >                                 dbl_mb(rhoge(1)),
1183     >                                 dbl_mb(rhoe(1)),
1184     >                                 dbl_mb(vctmp(1)),
1185     >                                 dbl_mb(vc(1)),
1186     >                                 dbl_mb(vcgg(1)),
1187     >                                 dbl_mb(vcge(1)),
1188     >                                 dbl_mb(vcee(1)))
1189         else if (ci_algorithm.eq.4) then
1190         call psi_4x4_virtual_gradient(dcpl_mb(psig_ptr),
1191     >                                 dbl_mb(psig_r(1)),
1192     >                                 dcpl_mb(psie_ptr),
1193     >                                 dbl_mb(psie_r(1)),
1194     >                                 eci,dcpl_mb(g(1)),
1195     >                                 dcpl_mb(t(1)),
1196     >                                 dcpl_mb(te(1)),
1197     >                                 dbl_mb(rhog(1)),
1198     >                                 dbl_mb(rhoge(1)),
1199     >                                 dbl_mb(rhoe(1)),
1200     >                                 dbl_mb(vctmp(1)),
1201     >                                 dbl_mb(vc(1)),
1202     >                                 dbl_mb(vcgg(1)),
1203     >                                 dbl_mb(vcge(1)),
1204     >                                 dbl_mb(vcee(1)))
1205         else if (ci_algorithm.eq.5) then
1206         call psi_2x2ne_virtual_gradient(dcpl_mb(psig_ptr),
1207     >                          dbl_mb(electron_psi_r_ptr()),
1208     >                          dcpl_mb(psig_homo(1)),
1209     >                          dbl_mb(psig_r(1)),
1210     >                          dcpl_mb(psie_ptr),
1211     >                          dbl_mb(psie_r(1)),
1212     >                          eci,dcpl_mb(g(1)),
1213     >                          dcpl_mb(t(1)),
1214     >                          dbl_mb(rhog(1)),
1215     >                          dbl_mb(vcee(1)),dbl_mb(vcge(1)),
1216     >                          dcpl_mb(hpsig(1)),dbl_mb(hpsig_r(1)))
1217         end if
1218
1219
1220         eold = eci
1221c         write(*,*) "it,eold=",it,eold
1222
1223c         if (oprint)
1224c     >   write(*,*) "orb,it,eci,eold,eci-eold=",ii,it,eci,eold,eci-eold
1225
1226         !*** project out to lower orbitals  ****
1227         call Pack_c_Copy(1,dcpl_mb(g(1)),dcpl_mb(r1(1)))
1228
1229*        **** preconditioning ****
1230         if (precondition) then
1231            pit = pit + 1
1232            call ke_Precondition(npack1,1,
1233     >                           dcpl_mb(psie_ptr),
1234     >                           dcpl_mb(r1(1)))
1235         end if
1236
1237         !call psi_project_out_virtual1(ii,dcpl_mb(r1(1)))
1238         call psi_project_out_virtual(ii,dcpl_mb(r1(1)))
1239
1240*        *** determine conjuagate direction ***
1241         call Pack_cc_dot(1,dcpl_mb(r1(1)),
1242     >                   dcpl_mb(r1(1)),
1243     >                   lmbda_r1)
1244
1245         call Pack_c_Copy(1,dcpl_mb(r1(1)),dcpl_mb(t(1)))
1246
1247         if (it.gt.1) then
1248            call Pack_cc_daxpy(1,(lmbda_r1/lmbda_r0),
1249     >                      dcpl_mb(t0(1)),
1250     >                      dcpl_mb(t(1)))
1251         end if
1252         lmbda_r0 = lmbda_r1
1253         oneloop = .true.
1254 3       call Pack_c_Copy(1,dcpl_mb(t(1)),dcpl_mb(t0(1)))
1255
1256*        *** normalize search direction, t ****
1257         call psi_project_out_virtual(ii,dcpl_mb(t(1)))
1258         call Pack_cc_dot(1,dcpl_mb(t(1)),dcpl_mb(t(1)),sigma)
1259         sigma = dsqrt(sigma)
1260         de0 = 1.0d0/sigma
1261         call Pack_c_SMul1(1,de0,dcpl_mb(t(1)))
1262
1263*        **** compute de0 = <t|g> ****
1264         call Pack_cc_dot(1,dcpl_mb(t(1)),
1265     >                      dcpl_mb(g(1)),
1266     >                      de0)
1267
1268*        *** bad direction ***
1269         if ((de0.lt.0.0d0).and.oneloop) then
1270           call Pack_c_Copy(1,dcpl_mb(g(1)),dcpl_mb(t(1)))
1271           oneloop = .false.
1272           go to 3
1273         end if
1274
1275         call psi_linesearch_CI_virtual(ci_algorithm,maxit_ls,
1276     >                                   ii,theta,eci,de0,dcpl_mb(t(1)),
1277     >                                  dcpl_mb(psig_ptr),
1278     >                                  dbl_mb(psig_r(1)),
1279     >                                  dbl_mb(psie_r(1)),
1280     >                                  dcpl_mb(r1(1)),dcpl_mb(te(1)),
1281     >                                  dbl_mb(rhog(1)),
1282     >                                  dbl_mb(rhoe(1)),
1283     >                                  dbl_mb(rhoge(1)),
1284     >                                  dbl_mb(vctmp(1)),
1285     >                                  dbl_mb(vc(1)),
1286     >                                  dbl_mb(vcgg(1)),
1287     >                                  dbl_mb(vcee(1)),
1288     >                                  dbl_mb(vcge(1)),
1289     >                                  psig_homo(1),
1290     >                                  electron_psi_r_ptr(),
1291     >                                  hpsig(1),hpsig_r(1),
1292     >                                  maxerror)
1293
1294         !precondition = (dabs(eci-eold).gt.(sp*maxerror))
1295
1296         if (oprint)
1297     >   write(*,*) "    orb,it,eci,eold,eci-eold=",
1298     >              ii,it,eci,eold,eci-eold,precondition
1299
1300         done = ((it.gt.maxiteration)
1301     >           .or.
1302     >           (dabs(eci-eold).lt.maxerror))
1303
1304         !write(*,*) "done=",done,it,maxiteration,dabs(eci-eold),maxerror
1305
1306         if (done) go to 4
1307      go to 2
1308
1309 4    e0old = e0
1310      call psi_get_gradient_virtual(ii,dcpl_mb(g(1)))
1311      call Pack_cc_dot(1,dcpl_mb(psie_ptr),dcpl_mb(g(1)),e0)
1312      e0 = e0
1313      error_out = dabs(e0-e0old)
1314
1315*     **** release stack memory ****
1316      value = .true.
1317      if (ci_algorithm.eq.5) then
1318         value = BA_pop_stack(hpsig_r(2))
1319         value = value.and.BA_pop_stack(hpsig(2))
1320         value = value.and.BA_pop_stack(psig_homo(2))
1321      end if
1322      value = value.and.BA_pop_stack(vcee(2))
1323      value = value.and.BA_pop_stack(vcge(2))
1324      value = value.and.BA_pop_stack(vcgg(2))
1325      value = value.and.BA_pop_stack(vc(2))
1326      value = value.and.BA_pop_stack(vctmp(2))
1327      value = value.and.BA_pop_stack(rhoge(2))
1328      value = value.and.BA_pop_stack(rhoe(2))
1329      value = value.and.BA_pop_stack(rhog(2))
1330      value = value.and.BA_pop_stack(psie_r(2))
1331      value = value.and.BA_pop_stack(psig_r(2))
1332      value = value.and.BA_pop_stack(te(2))
1333      value = value.and.BA_pop_stack(t(2))
1334      value = value.and.BA_pop_stack(g(2))
1335      value = value.and.BA_pop_stack(r1(2))
1336      value = value.and.BA_pop_stack(t0(2))
1337      if (.not. value) call errquit(
1338     >     'psi_2x2_update_virtual: popping stack memory',1, MA_ERR)
1339
1340      if (oprint) then
1341       write(luout,921) ii,e0,error_out,eci,dabs(eci-eold),it,pit,ep,sp
1342  921  format(5x,"orbital",I4," current e=",E14.6," error=",E9.3,
1343     >       " (eci=",E17.9,
1344     >       " error=",E9.3,")",
1345     >       " iterations",I4,"(",I4,
1346     >       " preconditioned, Ep,Sp=",F5.1,F7.1,")")
1347      end if
1348
1349
1350      return
1351      end
1352
1353
1354*     ***********************************
1355*     *				        *
1356*     *	     psi_linesearch_CI_virtual  *
1357*     *				        *
1358*     ***********************************
1359
1360*    This routine performs a linesearch on orbital ii, in the direction t.
1361* This routine is needed for a KS minimizer.
1362*  e0 = <orb|g>
1363*  de0 = 2*<t|g>
1364*
1365      subroutine psi_linesearch_CI_virtual(ci_algorithm,maxit,
1366     >                                     ii,theta,e0,de0,t,
1367     >                                     psig,psig_r,psie_r,
1368     >                                     h1psig,h1psie,
1369     >                                     rho,rhoge,rhoee,
1370     >                                     vctmp,vc,
1371     >                                     vcgg,vcee,vcge,
1372     >                                     psig_homo_ptr,psiall_r_ptr,
1373     >                                     hpsig_ptr,hpsig_r_ptr,
1374     >                                     maxerror)
1375      implicit none
1376      integer ci_algorithm,maxit
1377      integer ii
1378      real*8  theta
1379      real*8  e0,de0
1380      complex*16  t(*) !search direction
1381      complex*16  psig(*)
1382      real*8      psig_r(*)
1383      real*8      psie_r(*)
1384      complex*16  h1psig(*),h1psie(*)
1385      real*8      rho(*),rhoge(*),rhoee(*)
1386      real*8      vctmp(*),vc(*)
1387      real*8      vcgg(*),vcee(*),vcge(*)
1388      integer     psig_homo_ptr,psiall_r_ptr,hpsig_ptr,hpsig_r_ptr
1389      real*8      maxerror
1390
1391#include "errquit.fh"
1392#include "bafdecls.fh"
1393#include "psi.fh"
1394
1395*     **** local variables ****
1396      !integer maxit
1397      !parameter (maxit=15)
1398
1399      logical value,done,first
1400      integer orb(2),g(2),psi_ptr,it,iteration
1401      real*8 x,y,pi,dtheta_min,e1,de1,de1a,de1b
1402      real*8 theta0,theta1,theta2,dtheta,emin,theta_min,demin,tnorm
1403      real*8 ea,eb,ec,e2,de2
1404      real*8 emid,demid,thetatmp
1405
1406      psi_ptr=psi1_excited(1)
1407
1408      pi = 4.0d0*datan(1.0d0)
1409      dtheta_min = 0.01*theta
1410
1411*     **** allocate stack memory ****
1412      value = BA_push_get(mt_dcpl,npack1,'orb',
1413     >                       orb(2),orb(1))
1414      value = value.and.
1415     >        BA_push_get(mt_dcpl,npack1,'gg2',
1416     >                       g(2),g(1))
1417      if (.not. value) call errquit(
1418     >     'psi_linesearch_CI_virtual: out of stack memory',0, MA_ERR)
1419
1420
1421      call Pack_c_Copy(1,dcpl_mb(psi_ptr+(ii-1)*npack1),
1422     >                   dcpl_mb(orb(1)))
1423
1424      de0 = -de0
1425      call Pack_cc_dot(1,dcpl_mb(orb(1)),t,tnorm)
1426
1427      theta0 = 0.0d0
1428      iteration = 0
1429c      write(*,*)
1430c      write(*,*) "e0,de0, <psi|t> =",e0,de0,tnorm
1431c      write(*,*) "start line search:"
1432      emin = e0
1433      demin = de0
1434      dtheta = (2.0d0*pi/100.0d0)/99.0
1435      first = .true.
1436
1437c      if (ii.eq.2) then
1438c         write(*,*)
1439c         write(*,*) "e0,de0, <psi|t> =",e0,de0,tnorm
1440c         write(*,*) "start line search:"
1441c         dtheta = 0.05d0
1442c         thetatmp = theta
1443c         theta = 0.0d0
1444c      end if
1445
1446 10   iteration = iteration + 1
1447      x = cos(theta)
1448      y = sin(theta)
1449      call Pack_c_SMul(1,x,
1450     >                  dcpl_mb(orb(1)),
1451     >                  dcpl_mb(psi_ptr+(ii-1)*npack1))
1452       call Pack_cc_daxpy(1,y,
1453     >                   t,
1454     >                   dcpl_mb(psi_ptr+(ii-1)*npack1))
1455
1456*     *** determine e1,g1, and de1/dtheta ***
1457*     **** de/dtheta = -sin(theta)*<orb|g> + cos(theta)*<t|g> ****
1458      if (ci_algorithm.eq.2) then
1459      call psi_2x2_virtual_gradient(psig,psig_r,
1460     >                              dcpl_mb(psi_ptr+(ii-1)*npack1),
1461     >                              psie_r,
1462     >                              e1,dcpl_mb(g(1)),
1463     >                              h1psie,rho,vcee,vcge)
1464      else if (ci_algorithm.eq.3) then
1465      call psi_3x3_virtual_gradient(psig,psig_r,
1466     >                              dcpl_mb(psi_ptr+(ii-1)*npack1),
1467     >                              psie_r,
1468     >                              e1,dcpl_mb(g(1)),
1469     >                              h1psig,h1psie,
1470     >                              rho,rhoge,rhoee,
1471     >                              vctmp,vc,
1472     >                              vcgg,vcge,vcee)
1473      else if (ci_algorithm.eq.4) then
1474      call psi_4x4_virtual_gradient(psig,psig_r,
1475     >                              dcpl_mb(psi_ptr+(ii-1)*npack1),
1476     >                              psie_r,
1477     >                              e1,dcpl_mb(g(1)),
1478     >                              h1psig,h1psie,
1479     >                              rho,rhoge,rhoee,
1480     >                              vctmp,vc,
1481     >                              vcgg,vcge,vcee)
1482      else if (ci_algorithm.eq.5) then
1483      call psi_2x2ne_virtual_gradient(psig,dbl_mb(psiall_r_ptr),
1484     >                              dcpl_mb(psig_homo_ptr),psig_r,
1485     >                              dcpl_mb(psi_ptr+(ii-1)*npack1),
1486     >                              psie_r,
1487     >                              e1,dcpl_mb(g(1)),
1488     >                              h1psie,rho,vcee,vcge,
1489     >                              dcpl_mb(hpsig_ptr),
1490     >                              dbl_mb(hpsig_r_ptr))
1491      end if
1492      call Pack_cc_dot(1,dcpl_mb(orb(1)),dcpl_mb(g(1)),de1a)
1493      call Pack_cc_dot(1,t,              dcpl_mb(g(1)),de1b)
1494      de1 = y*de1a - x*de1b
1495c      write(*,*) "theta,e0,e1,de0,de1 = ",
1496c     >           ii,iteration,theta,e0,e1,de0,de1
1497
1498c      if (e1.lt.emin) then
1499c         emin  = e1
1500c         demin = de1
1501c         theta_min = theta
1502c      end if
1503
1504c      if ((first).and.(ii.eq.2)) then
1505c         write(*,*) "theta,e0,e1,de0,de1 = ",
1506c     >           ii,iteration,theta,e0,e1,de0,de1
1507c         theta = theta + dtheta
1508c         if (theta.gt.2.0d0*pi) then
1509c            first = .false.
1510c           dtheta = (2.0d0*pi/100.0d0)/99.0
1511c           theta = thetatmp
1512c
1513c         end if
1514c
1515c         go to 10
1516c      end if
1517
1518      if (e1.lt.emin) then
1519         emin  = e1
1520         demin = de1
1521         theta_min = theta
1522      end if
1523
1524      done = (de0.le.0.0d0).and.(de1.ge.0.0d0).or.(iteration.gt.maxit)
1525
1526      if (.not.done) then
1527         theta = theta + theta
1528         go to 10
1529      end if
1530
1531      theta0 = 0.0d0
1532      theta1 = theta
1533
1534      theta = theta1 - de0*(theta1-theta0)/(de1-de0)
1535
1536cccccc 20   theta = 0.5d0*(theta0+theta1)
1537 20   iteration = iteration + 1
1538      theta1 = theta
1539      x = cos(theta1)
1540      y = sin(theta1)
1541      call Pack_c_SMul(1,x,
1542     >                  dcpl_mb(orb(1)),
1543     >                  dcpl_mb(psi_ptr+(ii-1)*npack1))
1544       call Pack_cc_daxpy(1,y,
1545     >                   t,
1546     >                   dcpl_mb(psi_ptr+(ii-1)*npack1))
1547
1548*     *** determine e1,g1, and de1/dtheta ***
1549*     **** de/dtheta = -sin(theta)*<orb|g> + cos(theta)*<t|g> ****
1550      if (ci_algorithm.eq.2) then
1551      call psi_2x2_virtual_gradient(psig,psig_r,
1552     >                              dcpl_mb(psi_ptr+(ii-1)*npack1),
1553     >                              psie_r,
1554     >                              e1,dcpl_mb(g(1)),
1555     >                              h1psie,rho,vcee,vcge)
1556      else if (ci_algorithm.eq.3) then
1557      call psi_3x3_virtual_gradient(psig,psig_r,
1558     >                              dcpl_mb(psi_ptr+(ii-1)*npack1),
1559     >                              psie_r,
1560     >                              e1,dcpl_mb(g(1)),
1561     >                              h1psig,h1psie,
1562     >                              rho,rhoge,rhoee,
1563     >                              vctmp,vc,
1564     >                              vcgg,vcge,vcee)
1565      else if (ci_algorithm.eq.4) then
1566      call psi_4x4_virtual_gradient(psig,psig_r,
1567     >                              dcpl_mb(psi_ptr+(ii-1)*npack1),
1568     >                              psie_r,
1569     >                              e1,dcpl_mb(g(1)),
1570     >                              h1psig,h1psie,
1571     >                              rho,rhoge,rhoee,
1572     >                              vctmp,vc,
1573     >                              vcgg,vcge,vcee)
1574      else if (ci_algorithm.eq.5) then
1575      call psi_2x2ne_virtual_gradient(psig,dbl_mb(psiall_r_ptr),
1576     >                              dcpl_mb(psig_homo_ptr),psig_r,
1577     >                              dcpl_mb(psi_ptr+(ii-1)*npack1),
1578     >                              psie_r,
1579     >                              e1,dcpl_mb(g(1)),
1580     >                              h1psie,rho,vcee,vcge,
1581     >                              dcpl_mb(hpsig_ptr),
1582     >                              dbl_mb(hpsig_r_ptr))
1583      end if
1584      call Pack_cc_dot(1,dcpl_mb(orb(1)),dcpl_mb(g(1)),de1a)
1585      call Pack_cc_dot(1,t,              dcpl_mb(g(1)),de1b)
1586      de1 = y*de1a - x*de1b
1587      !dtheta =  -de0*(theta-theta0)/(demid-de0)
1588
1589      !write(*,*) "theta,e0,e1,(e1-e0)=",iteration,theta,e0,e1,e1-e0
1590      !write(*,*) "            de0,de1=",de0,de1
1591      if (e1.lt.emin) then
1592         emin  = e1
1593         demin = de1
1594         theta_min = theta1
1595      end if
1596
1597      if (e1.le.e0) then
1598         theta  = theta1 - de0*(theta1-theta0)/(de1-de0)
1599         e0     = e1
1600         de0    = de0
1601         theta0 = theta1
1602         done = .false.
1603      else if (e1.gt.e0) then
1604         theta = theta0 + 0.5d0*(theta1-theta0)
1605         done = .false.
1606      end if
1607
1608      done = (dabs(e1-e0).lt.maxerror).or.(iteration.gt.maxit)
1609
1610      if (.not.done) go to 20
1611
1612
1613
1614
1615c         theta = -pi/100.0 + (it-1)*(2.0d0*pi/100.0d0)/99.0
1616c
1617c*        **** orb2 = orb*cos(pi/300) + t*sin(pi/300) ****
1618c  10     x = cos(theta)
1619c         y = sin(theta)
1620c         call Pack_c_SMul(1,x,
1621c     >                  dcpl_mb(orb(1)),
1622c     >                  dcpl_mb(psi_ptr+(ii-1)*npack1))
1623c         call Pack_cc_daxpy(1,y,
1624c     >                   t,
1625c     >                   dcpl_mb(psi_ptr+(ii-1)*npack1))
1626c
1627c*        *** determine theta ***
1628c         call psi_2x2_virtual_gradient(psig,psig_r,
1629c     >                              dcpl_mb(psi_ptr+(ii-1)*npack1),
1630c     >                              psie_r,
1631c     >                              e1,dcpl_mb(g(1)),
1632c     >                              h1psie,rho,vcee,vcge)
1633c
1634c*        **** compute de/dtheta = -sin(theta)*<orb|g> + cos(theta)*<t|g> ****
1635c         call Pack_cc_dot(1,dcpl_mb(orb(1)),dcpl_mb(g(1)),de1a)
1636c         call Pack_cc_dot(1,t,              dcpl_mb(g(1)),de1b)
1637c         de1 = y*de1a - x*de1b
1638c
1639c         if (e1<emin) then
1640c            emin  = e1
1641c            demin = de1
1642c            theta_min = theta
1643c         end if
1644c
1645c
1646c         write(*,*) theta,e1,de1,(e1-ea)/(2.0d0*dtheta)
1647c         ea = eb
1648c         eb = e1
1649c      end do
1650c      theta = theta_min
1651      !write(*,*) ":end line search"
1652      !write(*,*)
1653      !write(*,*) "theta_min,emin,demin = ",theta_min,emin,demin
1654      !write(*,*)
1655
1656c*     **** compute de/dtheta = -sin(theta)*<orb|g> + cos(theta)*<t|g> ****
1657c      call Pack_cc_dot(1,dcpl_mb(orb(1)),dcpl_mb(g(1)),de1a)
1658c      call Pack_cc_dot(1,t,              dcpl_mb(g(1)),de1b)
1659c      de1 = y*de1a - x*de1b
1660c
1661c      write(*,*) "GERDD,ii,it,theta,e0,e1=",ii,it,theta,e0,e1,de0,de1
1662c     >
1663c
1664c      it = it + 1
1665c      if ((e1>e0).and.(it<5)) then
1666c         theta = 0.5d0*theta
1667c         goto 10
1668c      end if
1669c      dtheta = de0*(theta-0.0d0)/(de1-de0)
1670c      theta = theta + dtheta
1671c      write(*,*) "NEW theta = ",theta
1672
1673
1674
1675c      x = (e0 - e1 + 0.5d0*de0*sin(2*theta))
1676c     >    /(1.0d0-cos(2*theta))
1677c      theta = 0.5d0*datan(0.5d0*de0/x)
1678c      !write(*,*) "GERE,x,theta,e0,e1=",x,theta,e0,e1
1679
1680c*     **** orb2 = orb*cos(theta) + t*sin(theta) ****
1681c      x = cos(theta)
1682c      y = sin(theta)
1683c      call Pack_c_SMul(1,x,
1684c     >                  dcpl_mb(orb(1)),
1685c     >                  dcpl_mb(psi_ptr+(ii-1)*npack1))
1686c      call Pack_cc_daxpy(1,y,
1687c     >                   t,
1688c     >                   dcpl_mb(psi_ptr+(ii-1)*npack1))
1689
1690
1691*     **** release stack memory ****
1692      value =           BA_pop_stack(g(2))
1693      value = value.and.BA_pop_stack(orb(2))
1694      if (.not. value) call errquit(
1695     >     'psi_linesearch_CI_virtual: popping stack memory',1, MA_ERR)
1696
1697      return
1698      end
1699
1700c      subroutine psi_follow_geovirt(orb1,t,orb2)
1701c      implicit none
1702c      complex*16 orb(*),t(*),orb2(*)
1703c
1704c      integer i
1705c      real*8 theta x,y
1706c
1707c      do i=1,50
1708c         theta = (i-1)/50.0 * 2.0d0*datan(1.0d0)
1709c         call psi_linesearch_genpsi(theta,orb,t,orb2)
1710c
1711c         call psi_2x2_virtual_gradient(psig,psig_r,
1712c     >                              orb2,
1713c     >                              psie_r,
1714c     >                              e1,dcpl_mb(g(1)),
1715c     >                              h1psie,rho,vcee,vcge)
1716c      end do
1717
1718c      return
1719c      end
1720
1721
1722      subroutine psi_linesearch_genpsi(theta,orb,t,orb2)
1723      implicit none
1724      real*8 theta
1725      complex*16 orb(*),t(*),orb2(*)
1726      real*8 x,y
1727
1728      x = cos(theta)
1729      y = sin(theta)
1730      call Pack_c_SMul(1,x,orb,orb2)
1731      call Pack_cc_daxpy(1,y,t,orb2)
1732      return
1733      end
1734
1735*     ***********************************
1736*     *                                 *
1737*     *      psi_2x2_virtual_gradient   *
1738*     *                                 *
1739*     ***********************************
1740*
1741*    This routine calculates the 2x2 CI Energy and its gradient wrt to psie.
1742*
1743*    Entry - psig,psig_r
1744*            psie
1745*
1746*    Exit - E,dEpsie
1747*
1748*    Use - psie_r
1749*          h1psie,rho,vcee,vcge
1750*
1751      subroutine psi_2x2_virtual_gradient(psig,psig_r,psie,psie_r,
1752     >                                    E,dEdpsie,
1753     >                                    h1psie,rho,vcee,vcge)
1754      implicit none
1755      complex*16 psig(*)
1756      real*8     psig_r(*)
1757      complex*16 psie(*)
1758      real*8     psie_r(*)
1759      real*8     E
1760      complex*16 dEdpsie(*)
1761
1762      complex*16 h1psie(*)
1763      real*8     rho(*),vcee(*),vcge(*)
1764
1765#include "bafdecls.fh"
1766#include "psi.fh"
1767#include "errquit.fh"
1768#include "util.fh"
1769#include "stdio.fh"
1770
1771*     **** local variables ****
1772      logical oprint
1773      integer nx,ny,nz,n2ft3d
1774      real*8 Sgg,Sge,See,A,B,C,c1,c2
1775      real*8 H1gg,H1ge,H1ee
1776      real*8 H2gg,H2ge,H2ee
1777      real*8 Hgg,Hge,Hee
1778      real*8 ehartr,scal1,scal2,dv,Elow,Ehigh,lmbda
1779
1780*     **** external functions ****
1781      real*8   lattice_omega
1782      external lattice_omega
1783
1784      !write(*,*) "ENTER 2x2 GRADIENT"
1785      oprint = .false.
1786
1787      !*** generate psie_r ****
1788      call Pack_c_Copy(1,psie,psie_r)
1789      call Pack_c_unpack(1,psie_r)
1790      call D3dB_cr_fft3b(1,psie_r)
1791      call D3dB_r_Zero_Ends(1,psie_r)
1792
1793      call D3dB_n2ft3d(1,n2ft3d)
1794      call D3dB_nx(1,nx)
1795      call D3dB_ny(1,ny)
1796      call D3dB_nz(1,nz)
1797      scal1 = 1.0d0/dble(nx*ny*nz)
1798      scal2 = 1.0d0/lattice_omega()
1799      dv = scal1*lattice_omega()
1800
1801      !*** calculate Sgg, Sge, See ***
1802      call Pack_cc_dot(1,psig,psig,Sgg)
1803      call Pack_cc_dot(1,psig,psie,Sge)
1804      call Pack_cc_dot(1,psie,psie,See)
1805      !if (oprint) write(*,*) "Sgg=",Sgg," Sge=",Sge," See=",See
1806
1807*     **** apply H1 operator ****
1808      call  psi_H1psi(ispin,neq,npack1,n2ft3d,psig,psig_r,dEdpsie)
1809      call  psi_H1psi(ispin,neq,npack1,n2ft3d,psie,psie_r,h1psie)
1810
1811      call Pack_cc_dot(1,psig,dEdpsie,H1gg)
1812      call Pack_cc_dot(1,psig,h1psie,H1ge)
1813      call Pack_cc_dot(1,psie,h1psie,H1ee)
1814      H1gg = -H1gg
1815      H1ge = -H1ge
1816      H1ee = -H1ee
1817      if (ispin.eq.1) then
1818         H1gg = H1gg + H1gg
1819         H1ge = H1ge + H1ge
1820         H1ee = H1ee + H1ee
1821      end if
1822      H1ge = H1ge*Sge
1823      !if (oprint) write(*,*) "H1gg=",H1gg," H1ge=",H1ge," H1ee=",H1ee
1824
1825*     **** apply H2 operator  - note rho = 2*psi*psi ****
1826      call pspw_et_gen_rho(ispin,neq,n2ft3d,psig_r,rho)
1827      call coulomb2_v(rho,vcee)
1828      call D3dB_rr_dot(1,rho,vcee,ehartr)
1829      H2gg = 0.5d0*(0.5d0*ehartr*dv)
1830
1831      call pspw_et_gen_rho(ispin,neq,n2ft3d,psie_r,rho)
1832      call coulomb2_v(rho,vcee)
1833      call D3dB_rr_dot(1,rho,vcee,ehartr)
1834      H2ee = 0.5d0*(0.5d0*ehartr*dv)
1835
1836
1837      call pspw_et_gen_rho12(ispin,neq,n2ft3d,psig_r,psie_r,rho)
1838      call coulomb2_v(rho,vcge)
1839      call D3dB_rr_dot(1,rho,vcge,ehartr)
1840      H2ge = 0.5d0*(0.5d0*ehartr*dv)
1841
1842      !if (oprint) write(*,*) "H2gg=",H2gg," H2ge=",H2ge," H2ee=",H2ee
1843
1844*     **** generate and diagonalize 2x2 CI matrix ****
1845      Hgg = H1gg + H2gg
1846      Hge = H1ge + H2ge
1847      Hee = H1ee + H2ee
1848      A = 1.0d0
1849      B = -Hgg - Hee
1850      C = Hgg*Hee - Hge*Hge
1851      Elow  = (-B-dsqrt(B*B-4.0d0*A*C))/(2*A)
1852      Ehigh = (-B+dsqrt(B*B-4.0d0*A*C))/(2*A)
1853      lmbda = Hge/(Elow-Hee)
1854      c1 =  1.0d0/dsqrt(1.0d0+lmbda**2)
1855      !c2 =  dsqrt(1.0d0-c1*c1)
1856      c2 =  lmbda/dsqrt(1.0d0+lmbda**2)
1857      !if (oprint) write(*,*) "Elow=",Elow," Ehigh=",Ehigh
1858
1859*     **** generate dElow/dpsie ****
1860      call D3dB_rr_Mul2(1,psie_r,vcee)
1861      call D3dB_rr_Mul2(1,psig_r,vcge)
1862
1863      !call D3dB_r_SMul1(1,(-0.5d0*c2*c2),vcee)
1864      !call D3dB_rr_daxpy(1,(-0.5d0*c1*c2),vcge,vcee)
1865      call D3dB_r_SMul1(1,(-1.0d0*c2*c2),vcee)
1866      call D3dB_rr_daxpy(1,(-1.0d0*c1*c2),vcge,vcee)
1867
1868      call D3dB_rc_fft3f(1,vcee)
1869      call Pack_c_pack(1,vcee)
1870      call Pack_c_SMul(1,scal1,vcee,dEdpsie)
1871
1872      call Pack_cc_daxpy(1,2.0d0*c2*c2,h1psie,dEdpsie)
1873      !call Pack_cc_daxpy(1,2.0d0*c2*c2,h1psie,dEdpsie)
1874      E = Elow
1875      !if (oprint) write(*,*) "c1,c2=",c1,c2
1876
1877      return
1878      end
1879
1880*     ***********************************
1881*     *                                 *
1882*     *    psi_2x2ne_virtual_gradient   *
1883*     *                                 *
1884*     ***********************************
1885*
1886*    This routine calculates the 2x2 CI Energy and its gradient wrt to psie.
1887*
1888*    Entry - psig,psig_r
1889*            psie
1890*
1891*    Exit - E,dEpsie
1892*
1893*    Use - psie_r
1894*          h1psie,rho,vcee,vcge
1895*
1896      subroutine psi_2x2ne_virtual_gradient(psig,psig_r,
1897     >                                    psihomo,psihomo_r,
1898     >                                    psilumo,psilumo_r,
1899     >                                    E,dEdpsie,
1900     >                                    h1psie,rho,vcee,vcge,
1901     >                                    hpsig,hpsig_r)
1902      implicit none
1903      complex*16 psig(*)
1904      real*8     psig_r(*)
1905
1906      complex*16 psihomo(*)
1907      real*8     psihomo_r(*)
1908
1909      complex*16 psilumo(*)
1910      real*8     psilumo_r(*)
1911
1912      real*8     E
1913      complex*16 dEdpsie(*)
1914
1915      complex*16 h1psie(*)
1916      real*8     rho(*),vcee(*),vcge(*)
1917
1918      complex*16 hpsig(*)
1919      real*8     hpsig_r(*)
1920
1921#include "bafdecls.fh"
1922#include "psi.fh"
1923#include "errquit.fh"
1924#include "util.fh"
1925#include "stdio.fh"
1926
1927      real*8 Hgg
1928      common / CI_Hgg_common / Hgg
1929
1930*     **** local variables ****
1931      logical oprint
1932      integer neq1(2),i,q,ms,n1q(2),n2q(2)
1933      integer nx,ny,nz,n2ft3d
1934      real*8 Sgg,Sge,See,A,B,C,c1,c2
1935      real*8 H1gg,H1ge,H1ee
1936      real*8 H2gg,H2ge,H2ee
1937      real*8 Hge,Hee
1938      real*8 ehartr,scal1,scal2,dv,Elow,Ehigh,lmbda,sum1
1939      real*8 ehfx,phfx
1940
1941*     **** external functions ****
1942      real*8   lattice_omega
1943      external lattice_omega
1944
1945      !write(*,*) "ENTER 2x2 ne GRADIENT"
1946      oprint = .false.
1947
1948      !*** generate psilumo_r ****
1949      call Pack_c_Copy(1,psilumo,psilumo_r)
1950      call Pack_c_unpack(1,psilumo_r)
1951      call D3dB_cr_fft3b(1,psilumo_r)
1952      call D3dB_r_Zero_Ends(1,psilumo_r)
1953
1954      call D3dB_n2ft3d(1,n2ft3d)
1955      call D3dB_nx(1,nx)
1956      call D3dB_ny(1,ny)
1957      call D3dB_nz(1,nz)
1958      scal1 = 1.0d0/dble(nx*ny*nz)
1959      scal2 = 1.0d0/lattice_omega()
1960      dv = scal1*lattice_omega()
1961
1962      call Pack_c_Zero(1,dEdpsie)
1963
1964      !do i=1,5
1965      !   write(*,*) "i,psi_r,psi_homo=",psig_r(i),psihomo_r(i)
1966      !end do
1967
1968      !*** calculate Sgg, Sge, See ***
1969      Sgg = 1.0d0
1970      See = 1.0d0
1971      Sge = 0.0d0
1972      !if (oprint) write(*,*) "Sgg=",Sgg," Sge=",Sge," See=",See
1973
1974*     ***************************
1975*     **** one-electron part ****
1976*     ***************************
1977      n1q(1) = 1
1978      n2q(1) = neq(1)
1979      n1q(2) = neq(1)+1
1980      n2q(2) = neq(1)+neq(2)
1981
1982      !*** generate HF energy and gradient of state psig ****
1983      !*** generate Hgg if Hgg==0 ***
1984*     *** generate <g|H|g> ***
1985      if (dabs(Hgg).lt.1.0d-9) then
1986*        **** apply H1 operator ****
1987         call psi_H1psi(ispin,neq,npack1,n2ft3d,psig,psig_r,hpsig)
1988         H1gg=0.0d0
1989         do ms=1,ispin
1990         if (neq(ms).gt.0) then
1991            do q=n1q(ms),n2q(ms)
1992               call Pack_cc_idot(1,psig(1+npack1*(q-1)),
1993     >                            hpsig(1+npack1*(q-1)),sum1)
1994               H1gg = H1gg - sum1
1995            end do
1996         end if
1997         end do
1998         call Parallel_SumAll(H1gg)
1999         if (ispin.eq.1) H1gg = H1gg + H1gg
2000
2001         !write(*,*) "H1gg=",H1gg
2002
2003*        **** apply H2 operator  - note rho = 2*psi*psi ****
2004         H2gg=0.0d0
2005*        **** coulomb part ****
2006         call pspw_et_gen_rho(ispin,neq,n2ft3d,psig_r,rho)
2007         call coulomb2_v(rho,vcee)
2008         call D3dB_rr_dot(1,rho,vcee,ehartr)
2009         H2gg = (0.5d0*ehartr*dv)
2010         !write(*,*) "H2gg=",H2gg
2011
2012*        **** exchange part ****
2013         call D3dB_r_nZero(1,neq(1)+neq(2),hpsig_r)
2014         call pspw_potential_HFX(ispin,psig_r,hpsig_r)
2015         call pspw_energy_HFX(ispin,psig_r,ehfx,phfx)
2016         H2gg = H2gg + ehfx
2017         Hgg = H1gg + H2gg
2018         !write(*,*) "H2gg+ehfx=",H2gg
2019         !write(*,*) "Hgg=",Hgg
2020      end if
2021      !write(*,*) "Hgg=",Hgg
2022      !write(*,*)
2023
2024
2025
2026      !*** generate HF energy and gradient of state psie ****
2027*     *** generate <e|H|e> ***
2028      call Pack_c_Copy(1,psilumo,psig)
2029      call D3dB_c_Copy(1,psilumo_r,psig_r)
2030
2031*     **** apply H1 operator ****
2032      call psi_H1psi(ispin,neq,npack1,n2ft3d,psig,psig_r,hpsig)
2033      call Pack_c_Copy(1,hpsig,h1psie)
2034
2035      H1ee=0.0d0
2036      do ms=1,ispin
2037         if (neq(ms).gt.0) then
2038            do q=n1q(ms),n2q(ms)
2039               call Pack_cc_idot(1,psig(1+npack1*(q-1)),
2040     >                            hpsig(1+npack1*(q-1)),sum1)
2041               H1ee = H1ee - sum1
2042            end do
2043         end if
2044      end do
2045      call Parallel_SumAll(H1ee)
2046      if (ispin.eq.1) H1ee = H1ee + H1ee
2047      !write(*,*) "H1ee=",H1ee
2048
2049*     **** apply H2 operator  - note rho = 2*psi*psi ****
2050      H2ee=0.0d0
2051*     **** coulomb part ****
2052      call pspw_et_gen_rho(ispin,neq,n2ft3d,psig_r,rho)
2053      call coulomb2_v(rho,vcee)
2054      call D3dB_rr_dot(1,rho,vcee,ehartr)
2055      H2ee = (0.5d0*ehartr*dv)
2056      !write(*,*) "H2ee=",H2ee
2057
2058*     **** exchange part ****
2059      call D3dB_r_nZero(1,neq(1)+neq(2),hpsig_r)
2060      call pspw_potential_HFX(ispin,psig_r,hpsig_r)
2061      call pspw_energy_HFX(ispin,psig_r,ehfx,phfx)
2062      H2ee = H2ee + ehfx
2063      Hee = H1ee + H2ee
2064      !write(*,*) "H2ee+ehfx=",H2ee,ehfx
2065      !write(*,*) "Hee=",Hee
2066      !write(*,*)
2067
2068      call Pack_c_Copy(1,psihomo,psig)
2069      call D3dB_c_Copy(1,psihomo_r,psig_r)
2070
2071
2072*     *** generate <e|H|g> ***
2073      H1ge = 0.0d0
2074
2075      neq1(1) = 1
2076      neq1(2) = 0
2077      call pspw_et_gen_rho12(ispin,neq1,n2ft3d,psihomo_r,psilumo_r,rho)
2078      call coulomb2_v(rho,vcge)
2079      call D3dB_rr_dot(1,rho,vcge,ehartr)
2080      H2ge = 0.5d0*(0.5d0*ehartr*dv)
2081      Hge = H1ge + H2ge
2082      !write(*,*) "H1ge=",H1ge
2083      !write(*,*) "H2ge=",H2ge
2084      !write(*,*) "Hge=",Hge
2085      !write(*,*)
2086
2087*     **** generate and diagonalize 2x2 CI matrix ****
2088      A = 1.0d0
2089      B = -Hgg - Hee
2090      C = Hgg*Hee - Hge*Hge
2091      Elow  = (-B-dsqrt(B*B-4.0d0*A*C))/(2*A)
2092      Ehigh = (-B+dsqrt(B*B-4.0d0*A*C))/(2*A)
2093      lmbda = Hge/(Elow-Hee)
2094      c1 =  1.0d0/dsqrt(1.0d0+lmbda**2)
2095      !c2 =  dsqrt(1.0d0-c1*c1)
2096      c2 =  lmbda/dsqrt(1.0d0+lmbda**2)
2097      !write(*,*) "Elow=",Elow," Ehigh=",Ehigh
2098      !write(*,*)
2099
2100*     **** generate dElow/dpsie ****
2101      call D3dB_rr_Mul2(1,psilumo_r,vcee)
2102      call D3dB_rr_Sum2(1,hpsig_r,vcee)
2103
2104      call D3dB_rr_Mul2(1,psihomo_r,vcge)
2105
2106      call D3dB_r_SMul1(1,(-2.0d0*c2*c2),vcee)
2107      call D3dB_rr_daxpy(1,(-1.0d0*c1*c2),vcge,vcee)
2108      !original call D3dB_rr_daxpy(1,(-1.0d0*c1*c2),vcge,vcee)
2109
2110      call D3dB_rc_fft3f(1,vcee)
2111      call Pack_c_pack(1,vcee)
2112      call Pack_c_SMul(1,scal1,vcee,dEdpsie)
2113
2114      call Pack_cc_daxpy(1,2.0d0*c2*c2,h1psie,dEdpsie)
2115      E = Elow
2116      !if (oprint) write(*,*) "2x2ne c1,c2=",c1,c2
2117      write(*,*) "2x2ne c1,c2=",c1,c2,Elow
2118c
2119
2120      return
2121      end
2122
2123
2124
2125
2126*     ***********************************
2127*     *                                 *
2128*     *      psi_3x3_virtual_gradient   *
2129*     *                                 *
2130*     ***********************************
2131
2132*     This routine calculates the 3x3 CI Energy and its gradient wrt to psie.
2133*
2134*       g    e    m
2135*   g  Hgg  Hge  Hgm
2136*
2137*   e  Heg  Hee  Hem
2138*
2139*   m  Hmg  Hme  Hmm
2140
2141      subroutine psi_3x3_virtual_gradient(psig,psig_r,psie,psie_r,
2142     >                                    E,dEdpsie,
2143     >                                    h1psig,h1psie,
2144     >                                    rhogg,rhoge,rhoee,
2145     >                                    vctmp,vc,
2146     >                                    vcgg,vcge,vcee)
2147      implicit none
2148      complex*16 psig(*)
2149      real*8     psig_r(*)
2150      complex*16 psie(*)
2151      real*8     psie_r(*)
2152      real*8     E
2153      complex*16 dEdpsie(*)
2154
2155      complex*16 h1psie(*),h1psig(*)
2156      real*8     rhogg(*),rhoee(*),rhoge(*)
2157      real*8     vcgg(*),vcee(*),vcge(*)
2158      real*8     vctmp(*),vc(*)
2159
2160#include "bafdecls.fh"
2161#include "psi.fh"
2162#include "errquit.fh"
2163#include "util.fh"
2164#include "stdio.fh"
2165
2166*     **** local variables ****
2167      logical oprint
2168      integer nx,ny,nz,n2ft3d
2169      real*8 Sgg,Sge,See,c1,c2,c3,A,B,C
2170      real*8 H1gg,H1ge,H1gm,H1eg,H1ee,H1em,H1mg,H1me,H1mm
2171      real*8 H2gg,H2ge,H2gm,H2eg,H2ee,H2em,H2mg,H2me,H2mm
2172      real*8 Hgg,Hge,Hgm,Heg,Hee,Hem,Hmg,Hme,Hmm
2173      real*8 ehartr,scal1,scal2,dv,Elow,Ehigh,lmbda
2174      real*8 Hci(3,3),Eci(3),wkopt(1),work(999)
2175      !real*8 Hci(2,2),Eci(2),wkopt(1),work(999)
2176      integer lwork,INFO
2177
2178*     **** external functions ****
2179      real*8   lattice_omega
2180      external lattice_omega
2181
2182
2183      !write(*,*) "ENTER 3x3 GRADIENT"
2184      oprint = .false.
2185
2186      !*** generate psie_r ****
2187      call Pack_c_Copy(1,psie,psie_r)
2188      call Pack_c_unpack(1,psie_r)
2189      call D3dB_cr_fft3b(1,psie_r)
2190      call D3dB_r_Zero_Ends(1,psie_r)
2191
2192      call D3dB_n2ft3d(1,n2ft3d)
2193      call D3dB_nx(1,nx)
2194      call D3dB_ny(1,ny)
2195      call D3dB_nz(1,nz)
2196      scal1 = 1.0d0/dble(nx*ny*nz)
2197      scal2 = 1.0d0/lattice_omega()
2198      dv = scal1*lattice_omega()
2199
2200      !*** calculate Sgg, Sge, See ***
2201      call Pack_cc_dot(1,psig,psig,Sgg)
2202      call Pack_cc_dot(1,psig,psie,Sge)
2203      call Pack_cc_dot(1,psie,psie,See)
2204      !if (oprint) write(*,*) "Sgg=",Sgg," Sge=",Sge," See=",See
2205
2206
2207*     **** apply H1 operator ****
2208      call  psi_H1psi(ispin,neq,npack1,n2ft3d,psig,psig_r,h1psig)
2209      call  psi_H1psi(ispin,neq,npack1,n2ft3d,psie,psie_r,h1psie)
2210
2211      call Pack_cc_dot(1,psig,h1psig,H1gg)
2212      call Pack_cc_dot(1,psig,h1psie,H1ge)
2213      call Pack_cc_dot(1,psie,h1psie,H1ee)
2214      H1gg = -H1gg
2215      H1ge = -H1ge
2216      H1ee = -H1ee
2217      H1eg =  H1ge
2218      H1gm = dsqrt(2.0d0)*H1ge
2219      H1mm = H1gg + H1ee
2220      H1em = H1gm
2221      H1mg = H1gm
2222      H1me = H1em
2223      if (ispin.eq.1) then
2224         H1gg = H1gg + H1gg
2225         H1ge = H1ge + H1ge
2226         H1ee = H1ee + H1ee
2227         H1eg = H1ge
2228      end if
2229      H1ge = H1ge*Sge
2230      H1eg = H1eg*Sge
2231      !if (oprint) write(*,*) "H1gg=",H1gg," H1ge=",H1ge," H1ee=",H1ee
2232
2233*     **** apply H2 operator ****
2234
2235      call pspw_et_gen_rho(ispin,neq,n2ft3d,psig_r,rhogg)
2236      call coulomb2_v(rhogg,vcgg)
2237      call pspw_et_gen_rho12(ispin,neq,n2ft3d,psig_r,psie_r,rhoge)
2238      call coulomb2_v(rhoge,vcge)
2239      call pspw_et_gen_rho(ispin,neq,n2ft3d,psie_r,rhoee)
2240      call coulomb2_v(rhoee,vcee)
2241
2242      call D3dB_rr_dot(1,rhogg,vcgg,ehartr)
2243      H2gg = 0.5d0*(0.5d0*ehartr*dv)
2244
2245      call D3dB_rr_dot(1,rhoge,vcge,ehartr)
2246      H2ge = 0.5d0*(0.5d0*ehartr*dv)
2247
2248      call D3dB_rr_dot(1,rhogg,vcge,ehartr)
2249      H2gm = dsqrt(2.0d0)*0.5d0*(0.5d0*ehartr*dv)
2250
2251      call D3dB_rr_dot(1,rhoge,vcge,ehartr)
2252      H2eg = 0.5d0*(0.5d0*ehartr*dv)
2253
2254      call D3dB_rr_dot(1,rhoee,vcee,ehartr)
2255      H2ee = 0.5d0*(0.5d0*ehartr*dv)
2256
2257      call D3dB_rr_dot(1,rhoee,vcge,ehartr)
2258      H2em = dsqrt(2.0d0)*0.5d0*(0.5d0*ehartr*dv)
2259
2260      call D3dB_rr_dot(1,rhoge,vcgg,ehartr)
2261      H2mg = dsqrt(2.0d0)*0.5d0*(0.5d0*ehartr*dv)
2262
2263      call D3dB_rr_dot(1,rhoge,vcee,ehartr)
2264      H2me = dsqrt(2.0d0)*0.5d0*(0.5d0*ehartr*dv)
2265
2266      call D3dB_rr_dot(1,rhogg,vcee,ehartr)
2267      H2mm = 0.5d0*(0.5d0*ehartr*dv)
2268      H2mm = H2mm + H2ge
2269      !if (oprint) write(*,*) "H2gg=",H2gg," H2ge=",H2ge," H2ee=",H2ee
2270
2271
2272*     **** generate and diagonalize 2x2 CI matrix ****
2273      Hgg = H1gg + H2gg
2274      Hge = H1ge + H2ge
2275      Hgm = H1gm + H2gm
2276      Heg = H1eg + H2eg
2277      Hee = H1ee + H2ee
2278      Hem = H1em + H2em
2279      Hmg = H1mg + H2mg
2280      Hme = H1me + H2me
2281      Hmm = H1mm + H2mm
2282
2283      Hci(1,1) = Hgg
2284      Hci(1,2) = Hge
2285      Hci(1,3) = Hgm
2286      Hci(2,1) = Heg
2287      Hci(2,2) = Hee
2288      Hci(2,3) = Hem
2289      Hci(3,1) = Hmg
2290      Hci(3,2) = Hme
2291      Hci(3,3) = Hmm
2292
2293      call DSYEV('V','U',3,Hci,3,Eci,wkopt,-1,INFO )
2294      lwork = wkopt(1)
2295      call DSYEV('V','U',3,Hci,3,Eci,work,lwork,INFO )
2296
2297      c1 = Hci(1,1)
2298      c2 = Hci(2,1)
2299      c3 = Hci(3,1)
2300
2301*     **** generate dElow/dpsie ****
2302      call D3dB_r_Copy(1,vcee,vc)
2303      call D3dB_rr_Mul2(1,psie_r,vc)
2304      call D3dB_r_SMul1(1,(-2.0d0*c2*c2),vc)
2305
2306      call D3dB_r_Copy(1,vcge,vctmp)
2307      call D3dB_rr_Mul2(1,psig_r,vctmp)
2308      call D3dB_rr_daxpy(1,(-2.0d0*c1*c2),vctmp,vc)
2309
2310      call D3dB_r_Copy(1,vcee,vctmp)
2311      call D3dB_rr_Mul2(1,psig_r,vctmp)
2312      call D3dB_rr_daxpy(1,(-dsqrt(2.0d0)*c2*c3),vctmp,vc)
2313
2314      call D3dB_r_Copy(1,vcge,vctmp)
2315      call D3dB_rr_Mul2(1,psie_r,vctmp)
2316      call D3dB_rr_daxpy(1,(-dsqrt(2.0d0)*c2*c3),vctmp,vc)
2317
2318      call D3dB_r_Copy(1,vcgg,vctmp)
2319      call D3dB_rr_Mul2(1,psig_r,vctmp)
2320      call D3dB_rr_daxpy(1,(-dsqrt(2.0d0)*c3*c1),vctmp,vc)
2321
2322      call D3dB_r_Copy(1,vcge,vctmp)
2323      call D3dB_rr_Mul2(1,psie_r,vctmp)
2324      call D3dB_rr_daxpy(1,(-dsqrt(2.0d0)*c3*c2),vctmp,vc)
2325
2326      call D3dB_r_Copy(1,vcgg,vctmp)
2327      call D3dB_rr_Mul2(1,psie_r,vctmp)
2328      call D3dB_rr_daxpy(1,(-c3*c3),vctmp,vc)
2329
2330      call D3dB_r_Copy(1,vcge,vctmp)
2331      call D3dB_rr_Mul2(1,psig_r,vctmp)
2332      call D3dB_rr_daxpy(1,(-c3*c3),vctmp,vc)
2333
2334      call D3dB_rc_fft3f(1,vc)
2335      call Pack_c_pack(1,vc)
2336      call Pack_c_SMul(1,scal1,vc,dEdpsie)
2337
2338      call Pack_cc_daxpy(1,4.0d0*c2*c2,h1psie,dEdpsie)
2339      call Pack_cc_daxpy(1,2.0d0*dsqrt(2.0d0)*c2*c3,h1psig,dEdpsie)
2340      call Pack_cc_daxpy(1,2.0d0*dsqrt(2.0d0)*c3*c1,h1psig,dEdpsie)
2341      call Pack_cc_daxpy(1,2.0d0*c3*c3,h1psie,dEdpsie)
2342
2343      E = Eci(1)
2344      !if (oprint) write(*,*) "c1,c2=",c1,c2
2345
2346      return
2347      end
2348
2349
2350
2351*     ***********************************
2352*     *                                 *
2353*     *      psi_4x4_virtual_gradient   *
2354*     *                                 *
2355*     ***********************************
2356
2357*     This routine calculates the 4x4 CI Energy and its gradient wrt to psie.
2358*
2359*       g    e    a    b
2360*   g  Hgg  Hge  Hga  Hgb
2361*
2362*   e  Heg  Hee  Hea  Heb
2363*
2364*   a  Hag  Hae  Haa  Hab
2365*
2366*   b  Hbg  Hbe  Hba  Hbb
2367
2368      subroutine psi_4x4_virtual_gradient(psig,psig_r,psie,psie_r,
2369     >                                    E,dEdpsie,
2370     >                                    h1psig,h1psie,
2371     >                                    rhogg,rhoge,rhoee,
2372     >                                    vctmp,vc,
2373     >                                    vcgg,vcge,vcee)
2374      implicit none
2375      complex*16 psig(*)
2376      real*8     psig_r(*)
2377      complex*16 psie(*)
2378      real*8     psie_r(*)
2379      real*8     E
2380      complex*16 dEdpsie(*)
2381
2382      complex*16 h1psie(*),h1psig(*)
2383      real*8     rhogg(*),rhoee(*),rhoge(*)
2384      real*8     vcgg(*),vcee(*),vcge(*)
2385      real*8     vctmp(*),vc(*)
2386
2387#include "bafdecls.fh"
2388#include "psi.fh"
2389#include "errquit.fh"
2390#include "util.fh"
2391#include "stdio.fh"
2392
2393*     **** local variables ****
2394      integer nx,ny,nz,n2ft3d
2395      real*8 Sgg,Sge,See,c1,c2,c3,c4
2396      real*8 H1gg,H1ge,H1ga,H1gb,H1ee,H1ea,H1eb,H1aa,H1ab,H1bb
2397      real*8 H2gg,H2ge,H2ga,H2gb,H2ee,H2ea,H2eb,H2aa,H2ab,H2bb
2398      real*8 Hgg,Hge,Hga,Hgb,Hee,Hea,Heb,Haa,Hab,Hbb
2399      real*8 ehartr,scal1,scal2,dv,Elow,Ehigh,lmbda
2400      real*8 Hci(4,4),ciE(4),wkopt(1),work(999)
2401      integer lwork,INFO
2402
2403*     **** external functions ****
2404      real*8   lattice_omega
2405      external lattice_omega
2406
2407      !write(*,*) "ENTER 4x4 GRADIENT"
2408
2409      !*** generate psie_r ****
2410      call Pack_c_Copy(1,psie,psie_r)
2411      call Pack_c_unpack(1,psie_r)
2412      call D3dB_cr_fft3b(1,psie_r)
2413      call D3dB_r_Zero_Ends(1,psie_r)
2414
2415
2416      call D3dB_n2ft3d(1,n2ft3d)
2417      call D3dB_nx(1,nx)
2418      call D3dB_ny(1,ny)
2419      call D3dB_nz(1,nz)
2420      scal1 = 1.0d0/dble(nx*ny*nz)
2421      scal2 = 1.0d0/lattice_omega()
2422      dv = scal1*lattice_omega()
2423
2424      !*** calculate Sgg, Sge, See ***
2425      call Pack_cc_dot(1,psig,psig,Sgg)
2426      call Pack_cc_dot(1,psig,psie,Sge)
2427      call Pack_cc_dot(1,psie,psie,See)
2428      !if (oprint) write(*,*) "Sgg=",Sgg," Sge=",Sge," See=",See
2429
2430*     **** apply H1 operator ****
2431      call  psi_H1psi(ispin,neq,npack1,n2ft3d,psig,psig_r,h1psig)
2432      call  psi_H1psi(ispin,neq,npack1,n2ft3d,psie,psie_r,h1psie)
2433
2434      call Pack_cc_dot(1,psig,h1psig,H1gg)
2435      call Pack_cc_dot(1,psig,h1psie,H1ge)
2436      call Pack_cc_dot(1,psie,h1psie,H1ee)
2437      H1gg = -H1gg
2438      H1ge = -H1ge
2439      H1ee = -H1ee
2440      H1ga =  H1ge
2441      H1gb = -H1ge
2442      H1ea =  H1ge
2443      H1eb = -H1ge
2444      H1aa =  H1gg + H1ee
2445      H1ab =  H1gg + H1ee
2446      H1bb =  H1gg + H1ee
2447      if (ispin.eq.1) then
2448         H1gg =  H1gg + H1gg
2449         H1ge =  H1ge + H1ge
2450         H1ee =  H1ee + H1ee
2451      end if
2452      H1ge = H1ge*Sge
2453      H1ab = H1ab*Sge
2454
2455*     **** apply H2 operator ****
2456      call pspw_et_gen_rho(ispin,neq,n2ft3d,psig_r,rhogg)
2457      call coulomb2_v(rhogg,vcgg)
2458      call pspw_et_gen_rho(ispin,neq,n2ft3d,psie_r,rhoee)
2459      call coulomb2_v(rhoee,vcee)
2460      call pspw_et_gen_rho12(ispin,neq,n2ft3d,psig_r,psie_r,rhoge)
2461      call coulomb2_v(rhoge,vcge)
2462
2463      call D3dB_rr_dot(1,rhogg,vcgg,ehartr)
2464      H2gg = 0.5d0*(0.5d0*ehartr*dv)
2465
2466      call D3dB_rr_dot(1,rhoge,vcge,ehartr)
2467      H2ge = 0.5d0*(0.5d0*ehartr*dv)
2468
2469      call D3dB_rr_dot(1,rhoge,vcgg,ehartr)
2470      H2ga = 0.5d0*(0.5d0*ehartr*dv)
2471
2472      H2gb = -H2ga
2473
2474      call D3dB_rr_dot(1,rhoee,vcee,ehartr)
2475      H2ee = 0.5d0*(0.5d0*ehartr*dv)
2476
2477      call D3dB_rr_dot(1,rhoge,vcee,ehartr)
2478      H2ea = 0.5d0*(0.5d0*ehartr*dv)
2479
2480      H2eb = -H2ea
2481
2482      call D3dB_rr_dot(1,rhogg,vcee,ehartr)
2483      H2aa = 0.5d0*(0.5d0*ehartr*dv)
2484
2485      H2ab = -H2ge
2486
2487      H2bb = H2aa
2488
2489      !if (oprint) write(*,*) "H2gg=",H2gg," H2ge=",H2ge," H2ee=",H2ee
2490
2491*     **** generate and diagonalize 4x4 CI matrix ****
2492      Hgg = H1gg + H2gg
2493      Hge = H1ge + H2ge
2494      Hga = H1ga + H2ga
2495      Hgb = H1gb + H2gb
2496      Hee = H1ee + H2ee
2497      Hea = H1ea + H2ea
2498      Heb = H1eb + H2eb
2499      Haa = H1aa + H2aa
2500      Hab = H1ab + H2ab
2501      Hbb = H1bb + H2bb
2502
2503      Hci(1,1) = Hgg
2504      Hci(1,2) = Hge
2505      Hci(1,3) = Hga
2506      Hci(1,4) = Hgb
2507      Hci(2,1) = Hci(1,2)
2508      Hci(2,2) = Hee
2509      Hci(2,3) = Hea
2510      Hci(2,4) = Heb
2511      Hci(3,1) = Hci(1,3)
2512      Hci(3,2) = Hci(2,3)
2513      Hci(3,3) = Haa
2514      Hci(3,4) = Hab
2515      Hci(4,1) = Hci(1,4)
2516      Hci(4,2) = Hci(2,4)
2517      Hci(4,3) = Hci(3,4)
2518      Hci(4,4) = Hbb
2519
2520c     write(*,*) "Matrix"
2521c     write(*,*) Hci(1,1),Hci(1,2),Hci(1,3),Hci(1,4)
2522c     write(*,*) Hci(2,1),Hci(2,2),Hci(2,3),Hci(2,4)
2523c     write(*,*) Hci(3,1),Hci(3,2),Hci(3,3),Hci(3,4)
2524c     write(*,*) Hci(4,1),Hci(4,2),Hci(4,3),Hci(4,4)
2525c
2526      call DSYEV('V','U',4,Hci,4,ciE,wkopt,-1,INFO )
2527      lwork = wkopt(1)
2528      call DSYEV('V','U',4,Hci,4,ciE,work,lwork,INFO )
2529
2530c     write(*,*) "Eigen"
2531c     write(*,*) Hci(1,1),Hci(1,2),Hci(1,3),Hci(1,4)
2532c     write(*,*) Hci(2,1),Hci(2,2),Hci(2,3),Hci(2,4)
2533c     write(*,*) Hci(3,1),Hci(3,2),Hci(3,3),Hci(3,4)
2534c     write(*,*) Hci(4,1),Hci(4,2),Hci(4,3),Hci(4,4)
2535c     write(*,*) ciE(1),ciE(2),ciE(3),ciE(4)
2536
2537      c1 = Hci(1,1)
2538      c2 = Hci(2,1)
2539      c3 = Hci(3,1)
2540      c4 = Hci(4,1)
2541
2542*     **** generate dElow/dpsie ****
2543      call D3dB_r_Copy(1,vcee,vc)
2544      call D3dB_rr_Mul2(1,psie_r,vc)
2545      call D3dB_r_SMul1(1,(-2.0d0*c2*c2),vc)
2546
2547      call D3dB_r_Copy(1,vcge,vctmp)
2548      call D3dB_rr_Mul2(1,psig_r,vctmp)
2549      call D3dB_rr_daxpy(1,(-2.0d0*c2*c1),vctmp,vc)
2550
2551      call D3dB_r_Copy(1,vcge,vctmp)
2552      call D3dB_rr_Mul2(1,psie_r,vctmp)
2553      call D3dB_rr_daxpy(1,(-c2*c3),vctmp,vc)
2554
2555      call D3dB_r_Copy(1,vcee,vctmp)
2556      call D3dB_rr_Mul2(1,psig_r,vctmp)
2557      call D3dB_rr_daxpy(1,(-c2*c3),vctmp,vc)
2558
2559      call D3dB_r_Copy(1,vcge,vctmp)
2560      call D3dB_rr_Mul2(1,psie_r,vctmp)
2561      call D3dB_rr_daxpy(1,(c2*c4),vctmp,vc)
2562
2563      call D3dB_r_Copy(1,vcee,vctmp)
2564      call D3dB_rr_Mul2(1,psig_r,vctmp)
2565      call D3dB_rr_daxpy(1,(c2*c4),vctmp,vc)
2566
2567      call D3dB_r_Copy(1,vcgg,vctmp)
2568      call D3dB_rr_Mul2(1,psig_r,vctmp)
2569      call D3dB_rr_daxpy(1,(-c3*c1),vctmp,vc)
2570
2571      call D3dB_r_Copy(1,vcge,vctmp)
2572      call D3dB_rr_Mul2(1,psie_r,vctmp)
2573      call D3dB_rr_daxpy(1,(-c3*c2),vctmp,vc)
2574
2575      call D3dB_r_Copy(1,vcgg,vctmp)
2576      call D3dB_rr_Mul2(1,psie_r,vctmp)
2577      call D3dB_rr_daxpy(1,(-c3*c3),vctmp,vc)
2578
2579      call D3dB_r_Copy(1,vcge,vctmp)
2580      call D3dB_rr_Mul2(1,psig_r,vctmp)
2581      call D3dB_rr_daxpy(1,(c3*c4),vctmp,vc)
2582
2583      call D3dB_r_Copy(1,vcgg,vctmp)
2584      call D3dB_rr_Mul2(1,psig_r,vctmp)
2585      call D3dB_rr_daxpy(1,(c4*c1),vctmp,vc)
2586
2587      call D3dB_r_Copy(1,vcge,vctmp)
2588      call D3dB_rr_Mul2(1,psie_r,vctmp)
2589      call D3dB_rr_daxpy(1,(c4*c2),vctmp,vc)
2590
2591      call D3dB_r_Copy(1,vcge,vctmp)
2592      call D3dB_rr_Mul2(1,psig_r,vctmp)
2593      call D3dB_rr_daxpy(1,(c4*c3),vctmp,vc)
2594
2595      call D3dB_r_Copy(1,vcgg,vctmp)
2596      call D3dB_rr_Mul2(1,psie_r,vctmp)
2597      call D3dB_rr_daxpy(1,(-c4*c4),vctmp,vc)
2598
2599      call D3dB_rc_fft3f(1,vc)
2600      call Pack_c_pack(1,vc)
2601      call Pack_c_SMul(1,scal1,vc,dEdpsie)
2602
2603      call Pack_cc_daxpy(1, 4.0d0*c2*c2,h1psie,dEdpsie)
2604      call Pack_cc_daxpy(1, 2.0d0*c2*c3,h1psig,dEdpsie)
2605      call Pack_cc_daxpy(1,-2.0d0*c2*c4,h1psig,dEdpsie)
2606      call Pack_cc_daxpy(1, 2.0d0*c3*c1,h1psig,dEdpsie)
2607      call Pack_cc_daxpy(1, 2.0d0*c3*c3,h1psie,dEdpsie)
2608      call Pack_cc_daxpy(1,-2.0d0*c4*c1,h1psig,dEdpsie)
2609      call Pack_cc_daxpy(1, 2.0d0*c4*c4,h1psie,dEdpsie)
2610
2611      E = ciE(1)
2612
2613      return
2614      end
2615
2616
2617
2618
2619
2620*     ****************************************************
2621*     *                                                  *
2622*     *                psi_H1psi                         *
2623*     *                                                  *
2624*     ****************************************************
2625      subroutine psi_H1psi(ispin,neq,npack1,n2ft3d,
2626     >                        psi,psi_r,
2627     >                        H1psi)
2628      implicit none
2629      integer ispin,neq(2),npack1,n2ft3d
2630      complex*16 psi(npack1,*)
2631      real*8     psi_r(n2ft3d,*)
2632      real*8     H1psi(npack1,*)
2633
2634#include "bafdecls.fh"
2635#include "errquit.fh"
2636
2637*     **** local variables ****
2638      logical value,aperiodic,field_exist
2639      integer q,n,m,ms,nfft3d,nemaxq,npack0,n1q(2),n2q(2)
2640      integer Hpsi2(2),vlr_l(2),r_grid(2),v_field(2),vl(2)
2641      real*8 tmp1(2),tmp2(2),sum1,sum2
2642
2643*     **** external functions ****
2644      logical  pspw_charge_found,pspw_Efield_found
2645      external pspw_charge_found,pspw_Efield_found
2646      integer  control_version
2647      external control_version
2648
2649      aperiodic   = (control_version().eq.4)
2650      field_exist = pspw_charge_found().or.pspw_Efield_found()
2651      nemaxq = neq(1)+neq(2)
2652
2653      call Pack_npack(0,npack0)
2654      call D3dB_nfft3d(1,nfft3d)
2655
2656*     **** allocate memory ****
2657      value = .true.
2658      if (aperiodic) then
2659       value = value.and.
2660     >        BA_push_get(mt_dbl,(n2ft3d),'vlr_l',vlr_l(2),vlr_l(1))
2661      end if
2662      if (field_exist.or.aperiodic) then
2663         value = value.and.
2664     >            BA_push_get(mt_dbl,(3*n2ft3d),'r_grid',
2665     >                        r_grid(2),r_grid(1))
2666         value = value.and.
2667     >           BA_push_get(mt_dbl,(n2ft3d),'v_field',
2668     >                       v_field(2),v_field(1))
2669      end if
2670      value = value.and.
2671     >        BA_push_get(mt_dcpl,(npack0),'vloc',vl(2),vl(1))
2672      if (.not. value)
2673     >   call errquit('pspw_et_sub1:out of stack',0,MA_ERR)
2674
2675
2676*     **** generate r_grid ****
2677      if (aperiodic.or.field_exist)
2678     >   call lattice_r_grid(dbl_mb(r_grid(1)))
2679
2680*     **** generate local pseudopotential  ****
2681      call v_local(dcpl_mb(vl(1)),.false.,tmp1,tmp2)
2682
2683*     *** long-range psp for charge systems ***
2684      if (control_version().eq.4) then
2685         call v_lr_local(dbl_mb(r_grid(1)),dbl_mb(vlr_l(1)))
2686      end if
2687
2688*     ***** generate other real-space fields ****
2689      if (field_exist) then
2690         !call dcopy(n2ft3d,0.0d0,0,dbl_mb(v_field(1)),1)
2691         call Parallel_shared_vector_zero(.true.,n2ft3d,
2692     >                                    dbl_mb(v_field(1)))
2693         call pspw_charge_Generate_V(n2ft3d,
2694     >                               dbl_mb(r_grid(1)),
2695     >                               dbl_mb(v_field(1)))
2696         call pspw_Efield_Generate_V(n2ft3d,
2697     >                               dbl_mb(r_grid(1)),
2698     >                               dbl_mb(v_field(1)))
2699      end if
2700
2701
2702*     **** get H1psi  ****
2703      !call dcopy(2*nemaxq*npack1,0.0d0,0,H1psi,1)
2704      call Parallel_shared_vector_zero(.true.,2*nemaxq*npack1,H1psi)
2705      if (aperiodic) then
2706         call psi_H1v4(ispin,neq,psi,psi_r,
2707     >             dcpl_mb(vl(1)),dbl_mb(vlr_l(1)),
2708     >             dbl_mb(v_field(1)),field_exist,
2709     >             H1psi)
2710      else
2711         call psi_H1(ispin,neq,psi,psi_r,
2712     >             dcpl_mb(vl(1)),
2713     >             dbl_mb(v_field(1)),field_exist,
2714     >             H1psi)
2715      end if
2716
2717
2718*     **** deallocate memory ****
2719      value = BA_pop_stack(vl(2))
2720      if (field_exist.or.aperiodic) then
2721         value = value.and.BA_pop_stack(v_field(2))
2722         value = value.and.BA_pop_stack(r_grid(2))
2723      end if
2724      if (aperiodic) then
2725         value = value.and.BA_pop_stack(vlr_l(2))
2726      end if
2727      if (.not.value)
2728     >   call errquit('psi_H1psi:pop stack',0,MA_ERR)
2729
2730      return
2731      end
2732
2733
2734
2735
2736
2737
2738
2739************************ virtural orbital Part ************************
2740*     ***********************************
2741*     *				        *
2742*     *	     psi_minimize_virtual       *
2743*     *				        *
2744*     ***********************************
2745      subroutine psi_minimize_virtual()
2746      implicit none
2747
2748#include "bafdecls.fh"
2749#include "psi.fh"
2750
2751      !*** local variables ***
2752      integer maxit_orb
2753      integer ii,l,l2
2754      real*8  sum,maxerror,error_out,e0
2755
2756      !*** external functions ***
2757      real*8   control_tole
2758      external control_tole
2759
2760      !call psi_gen_density_potentials(1)
2761      maxit_orb=120
2762      maxerror = control_tole()
2763
2764      do ii=1,(ne_excited(1)+ne_excited(2))
2765         l2= 0
2766
2767         !*** orthogonalize to lower orbitals  ****
2768 2       l2 = l2 + 1
2769         call psi_project_out_virtual1(
2770     >           ii,
2771     >           dcpl_mb(psi1_excited(1)+(ii-1)*npack1))
2772
2773         !*** normalize ****
2774         call Pack_cc_dot(1,
2775     >            dcpl_mb(psi1_excited(1) +(ii-1)*npack1),
2776     >            dcpl_mb(psi1_excited(1) +(ii-1)*npack1),
2777     >            sum)
2778         sum = 1.0d0/dsqrt(sum)
2779         call Pack_c_SMul1(1,sum,
2780     >            dcpl_mb(psi1_excited(1) +(ii-1)*npack1))
2781
2782         !*** minimize orbital ****
2783          l = 0
2784 3        call psi_KS_update_virtual(maxit_orb,
2785     >                               maxerror,
2786     >                               0.001d0,ii,error_out,e0)
2787          l  = l+1
2788          if ((error_out.gt.maxerror).and.(l.le.(1+(l2-1)*3))) go to 3
2789          if (((error_out.gt.maxerror).or.(e0.gt.4.0d0))
2790     >        .and.(l2.le.1)) then
2791            call Pack_c_Zero(1,
2792     >               dcpl_mb(psi1_excited(1) +(ii-1)*npack1))
2793            call Pack_c_setzero(1,1.0d0,
2794     >               dcpl_mb(psi1_excited(1) +(ii-1)*npack1))
2795             go to 2
2796          end if
2797
2798          dbl_mb(eig_excited(1)+ii-1) = e0
2799
2800      end do
2801      call psi_sort_virtual()
2802
2803
2804      return
2805      end
2806
2807      subroutine psi_sort_virtual()
2808      implicit none
2809
2810#include "bafdecls.fh"
2811#include "psi.fh"
2812#include "errquit.fh"
2813
2814      logical value
2815      integer i,j,ii,jj,ms
2816      integer r1(2)
2817      real*8  ei,ej
2818
2819      value = BA_push_get(mt_dcpl,npack1,'r1',r1(2),r1(1))
2820      if (.not. value) call errquit(
2821     >     'psi_sort_virtual: out of stack memory',0, MA_ERR)
2822
2823      do ms=1,ispin
2824
2825        !*** Bubble sort ***
2826        do ii=1,ne_excited(ms)
2827         do jj=ii+1,ne_excited(ms)
2828           i = ii + (ms-1)*ne_excited(1)
2829           j = jj + (ms-1)*ne_excited(1)
2830           ei = dbl_mb(eig_excited(1)+i-1)
2831           ej = dbl_mb(eig_excited(1)+j-1)
2832
2833           !*** swap ***
2834           if (ej.lt.ei) then
2835             dbl_mb(eig_excited(1)+i-1) = ej
2836             dbl_mb(eig_excited(1)+j-1) = ei
2837             call Pack_c_Copy(1,dcpl_mb(psi1_excited(1)+(i-1)*npack1),
2838     >                          dcpl_mb(r1(1)))
2839             call Pack_c_Copy(1,dcpl_mb(psi1_excited(1)+(j-1)*npack1),
2840     >                          dcpl_mb(psi1_excited(1)+(i-1)*npack1))
2841             call Pack_c_Copy(1,dcpl_mb(r1(1)),
2842     >                          dcpl_mb(psi1_excited(1)+(j-1)*npack1))
2843           end if
2844
2845         end do
2846        end do
2847
2848      end do
2849
2850      value = BA_pop_stack(r1(2))
2851      if (.not. value) call errquit(
2852     >     'psi_sort_virtual: popping stack memory',1, MA_ERR)
2853      return
2854      end
2855
2856*     ***********************************
2857*     *				        *
2858*     *	     psi_KS_update_virtual      *
2859*     *				        *
2860*     ***********************************
2861
2862*    This routine performs a KS update on virtual ii
2863*
2864      subroutine psi_KS_update_virtual(maxiteration,
2865     >                             maxerror,perror,ii,
2866     >                             error_out,e0)
2867      implicit none
2868      integer maxiteration
2869      real*8  maxerror,perror
2870      integer ii
2871      real*8 error_out
2872      real*8 e0
2873
2874#include "bafdecls.fh"
2875#include "psi.fh"
2876#include "errquit.fh"
2877#include "util.fh"
2878#include "stdio.fh"
2879
2880*     **** local variables ****
2881      integer MASTER,taskid
2882      parameter (MASTER=0)
2883
2884      logical value,done,oneloop,precondition,oprint
2885      integer it,pit
2886      real*8 eold,percent_error,error0,de0,lmbda_r0,lmbda_r1
2887      real*8 theta,ep,sp
2888      integer r1(2),t0(2),t(2),g(2)
2889      integer psi_ptr
2890
2891      logical  control_print
2892      external control_print
2893      real*8   control_Ep,control_Sp
2894      external control_Ep,control_Sp
2895
2896      psi_ptr=psi1_excited(1)
2897
2898      call Parallel_taskid(taskid)
2899      oprint= ((taskid.eq.MASTER).and.control_print(print_medium))
2900
2901
2902      value = BA_push_get(mt_dcpl,npack1,'t0',t0(2),t0(1))
2903      value = value.and.
2904     >        BA_push_get(mt_dcpl,npack1,'r1',r1(2),r1(1))
2905      value = value.and.
2906     >        BA_push_get(mt_dcpl,npack1,'g',g(2),g(1))
2907      value = value.and.
2908     >        BA_push_get(mt_dcpl,npack1,'t',t(2),t(1))
2909      if (.not. value) call errquit(
2910     >     'psi_KS_update_virtual: out of stack memory',0, MA_ERR)
2911
2912      ep = control_Ep()
2913      sp = control_Sp()
2914      precondition = .true.
2915      done = .false.
2916      error0 = 0.0d0
2917      e0 = 0.0d0
2918      theta = -3.14159d0/600.0d0
2919      lmbda_r0 = 1.0d0
2920      it = 0
2921      pit = 0
2922 2    continue
2923
2924         it = it + 1
2925         eold = e0
2926
2927*        *** calculate residual (steepest descent) direction for a single band ***
2928         call psi_get_gradient_virtual(ii,dcpl_mb(g(1)))
2929         call Pack_cc_dot(1,dcpl_mb(psi_ptr+(ii-1)*npack1),
2930     >                   dcpl_mb(g(1)),
2931     >                    e0)
2932
2933         e0 = -e0
2934
2935         !if (it.eq.1) then
2936         !  percent_error = 1.0d0
2937         !else if (it.eq.2) then
2938         !  error0 = dabs(e0-eold)
2939         !  percent_error = 1.0d0
2940         !else
2941         percent_error=0.0d0
2942         if(error0.ne.0.0d0)
2943     A      percent_error = dabs(e0-eold)/error0
2944         !end if
2945
2946         precondition = (dabs(e0-eold).gt.(sp*maxerror))
2947
2948         done = ((it.gt.maxiteration)
2949     >           .or.
2950     >           (dabs(e0-eold).lt.maxerror))
2951
2952         if (done) go to 4
2953
2954
2955         call Pack_c_Copy(1,dcpl_mb(g(1)),dcpl_mb(r1(1)))
2956         call Pack_cc_daxpy(1,(e0),
2957     >                 dcpl_mb(psi_ptr+(ii-1)*npack1),
2958     >                 dcpl_mb(r1(1)))
2959
2960*        **** preconditioning ****
2961         if (precondition) then
2962            pit = pit + 1
2963            call ke_Precondition(npack1,1,
2964     >                           dcpl_mb(psi_ptr+(ii-1)*npack1),
2965     >                           dcpl_mb(r1(1)))
2966         end if
2967
2968         !*** determine conjuagate direction ***
2969         call Pack_cc_dot(1,dcpl_mb(r1(1)),
2970     >                   dcpl_mb(r1(1)),
2971     >                   lmbda_r1)
2972         call Pack_c_Copy(1,dcpl_mb(r1(1)),dcpl_mb(t(1)))
2973
2974         if (it.gt.1) then
2975         call Pack_cc_daxpy(1,(lmbda_r1/lmbda_r0),
2976     >                   dcpl_mb(t0(1)),
2977     >                   dcpl_mb(t(1)))
2978         end if
2979         lmbda_r0 = lmbda_r1
2980         oneloop = .true.
2981 3       call Pack_c_Copy(1,dcpl_mb(t(1)),dcpl_mb(t0(1)))
2982
2983
2984*        *** normalize search direction, t ****
2985         call psi_project_out_virtual(ii,dcpl_mb(t(1)))
2986         call Pack_cc_dot(1,dcpl_mb(t(1)),
2987     >                   dcpl_mb(t(1)),
2988     >                   de0)
2989         de0 = 1.0d0/dsqrt(de0)
2990c         call Pack_c_SMul(1,de0,dcpl_mb(t(1)),dcpl_mb(t(1)))
2991         call Pack_c_SMul1(1,de0,dcpl_mb(t(1)))
2992         call Pack_cc_dot(1,dcpl_mb(t(1)),
2993     >                   dcpl_mb(g(1)),
2994     >                   de0)
2995
2996*        *** bad direction ***
2997         if ((de0.lt.0.0d0).and.oneloop) then
2998           call Pack_c_Copy(1,dcpl_mb(g(1)),dcpl_mb(t(1)))
2999           oneloop = .false.
3000           go to 3
3001         end if
3002
3003         de0 = -2.0d0*de0
3004         call psi_linesearch_virtual(ii,
3005     >                               theta,e0,de0,dcpl_mb(t(1)))
3006
3007      go to 2
3008
3009
3010*     **** release stack memory ****
3011 4    value =           BA_pop_stack(t(2))
3012      value = value.and.BA_pop_stack(g(2))
3013      value = value.and.BA_pop_stack(r1(2))
3014      value = value.and.BA_pop_stack(t0(2))
3015      if (.not. value) call errquit(
3016     >     'psi_KS_update_virtual: popping stack memory',1, MA_ERR)
3017
3018      if (oprint) then
3019         write(luout,921) ii,-e0,dabs(e0-eold),it,pit,ep,sp
3020  921 format(5x,"orbital",I4," current e=",E10.3,
3021     >       " (error=",E9.3,")",
3022     >       " iterations",I4,"(",I4,
3023     >       " preconditioned, Ep,Sp=",F5.1,F7.1,")")
3024      end if
3025
3026      error_out = dabs(e0-eold)
3027      e0        = -e0
3028      return
3029      end
3030
3031*     ***********************************
3032*     *				        *
3033*     *	     psi_linesearch_virtual     *
3034*     *				        *
3035*     ***********************************
3036
3037*    This routine performs a linesearch on orbital ii, in the direction t.
3038* This routine is needed for a KS minimizer.
3039*  e0 = <orb|g>
3040*  de0 = 2*<t|g>
3041*
3042      subroutine psi_linesearch_virtual(ii,theta,e0,de0,t)
3043      implicit none
3044      integer ii
3045      real*8  theta
3046      real*8  e0,de0
3047      complex*16 t(*) !search direction
3048#include "errquit.fh"
3049#include "bafdecls.fh"
3050#include "psi.fh"
3051
3052*     **** local variables ****
3053      logical value
3054      integer orb(2),g(2),psi_ptr
3055      real*8 x,y,pi,dtheta_min,e1
3056
3057      psi_ptr=psi1_excited(1)
3058
3059      pi = 4.0d0*datan(1.0d0)
3060      !dtheta = pi/300.0d0
3061      dtheta_min = 0.01*theta
3062
3063*     **** allocate stack memory ****
3064      value = BA_push_get(mt_dcpl,npack1,'orb',
3065     >                       orb(2),orb(1))
3066      value = value.and.
3067     >        BA_push_get(mt_dcpl,npack1,'g',
3068     >                       g(2),g(1))
3069      if (.not. value) call errquit(
3070     >     'psi_linesearch_virtual: out of stack memory',0, MA_ERR)
3071
3072
3073      call Pack_c_Copy(1,dcpl_mb(psi_ptr+(ii-1)*npack1),
3074     >                   dcpl_mb(orb(1)))
3075
3076*     **** orb2 = orb*cos(pi/300) + t*sin(pi/300) ****
3077  10  x = cos(theta)
3078      y = sin(theta)
3079      call Pack_c_SMul(1,x,
3080     >                  dcpl_mb(orb(1)),
3081     >                  dcpl_mb(psi_ptr+(ii-1)*npack1))
3082      call Pack_cc_daxpy(1,y,
3083     >                   t,
3084     >                   dcpl_mb(psi_ptr+(ii-1)*npack1))
3085
3086*     *** determine theta ***
3087      call psi_get_gradient_virtual(ii,dcpl_mb(g(1)))
3088      call Pack_cc_dot(1,dcpl_mb(psi_ptr+(ii-1)*npack1),
3089     >                   dcpl_mb(g(1)),
3090     >                   e1)
3091      e1 = -e1
3092
3093
3094      !if (((-e1).gt.(-e0)).and.(theta.gt.dtheta_min)) then
3095      !   theta = 0.5d0*theta
3096      !   go to 10
3097      !end if
3098
3099      x = (e0 - e1 + 0.5d0*de0*sin(2*theta))
3100     >    /(1.0d0-cos(2*theta))
3101      theta = 0.5d0*datan(0.5d0*de0/x)
3102
3103
3104
3105*     **** orb2 = orb*cos(theta) + t*sin(theta) ****
3106      x = cos(theta)
3107      y = sin(theta)
3108      call Pack_c_SMul(1,x,
3109     >                  dcpl_mb(orb(1)),
3110     >                  dcpl_mb(psi_ptr+(ii-1)*npack1))
3111      call Pack_cc_daxpy(1,y,
3112     >                   t,
3113     >                   dcpl_mb(psi_ptr+(ii-1)*npack1))
3114
3115
3116*     **** release stack memory ****
3117      value =           BA_pop_stack(g(2))
3118      value = value.and.BA_pop_stack(orb(2))
3119      if (.not. value) call errquit(
3120     >     'psi_linesearch_virtual: popping stack memory',1, MA_ERR)
3121
3122      return
3123      end
3124
3125
3126*     ***********************************
3127*     *				        *
3128*     *	     psi_get_gradient_virtual	*
3129*     *				        *
3130*     ***********************************
3131
3132*    This routine returns the Hpsi(i).
3133* This routine is needed for a KS minimizer.
3134*
3135      subroutine psi_get_gradient_virtual(ii,Horb)
3136      implicit none
3137      integer ii
3138      complex*16 Horb(*)
3139
3140#include "bafdecls.fh"
3141#include "psi.fh"
3142
3143*     **** local variables ****
3144      integer psi_ptr,ms
3145
3146      psi_ptr=psi1_excited(1)+(ii-1)*npack1
3147
3148      if (ii.le.ne_excited(1)) then
3149         ms = 1
3150      else
3151         ms = 2
3152      end if
3153
3154      call electron_get_gradient_virtual(ms,dcpl_mb(psi_ptr),Horb)
3155
3156      return
3157      end
3158
3159*     *******************************************
3160*     *				                *
3161*     *	         psi_project_out_virtual1        *
3162*     *				                *
3163*     *******************************************
3164*
3165*    This routine projects out non-orthogonal components of Horb.
3166* This routine is needed for a KS minimizer.
3167*
3168      subroutine psi_project_out_virtual1(ii,Horb)
3169      integer ii
3170      complex*16 Horb(*)
3171
3172#include "bafdecls.fh"
3173#include "errquit.fh"
3174#include "psi.fh"
3175
3176      integer ms,i,jj,kk,shift,shifte
3177      integer etmp(2)
3178      real*8  sum
3179
3180*     **** spin up orbital ****
3181      if (ii.le.ne_excited(1)) then
3182
3183         shift  = 0
3184         shifte = 0
3185         ms     = 1
3186         kk     = ii
3187*     **** spin down orbital ****
3188      else
3189         shift  = neq(1)*npack1
3190         shifte = ne_excited(1)*npack1
3191         ms     = 2
3192         kk     = ii-ne_excited(1)
3193      end if
3194
3195      !**** project out filled orbitals ****
3196      if (neq(ms).eq.ne(ms)) then
3197
3198         do i=1,ne(ms)
3199           call Pack_cc_dot(1,
3200     >            dcpl_mb(psi1(1) +(i-1)*npack1+shift),
3201     >            Horb,
3202     >            sum)
3203           call Pack_cc_daxpy(1,(-sum),
3204     >            dcpl_mb(psi1(1) +(i-1)*npack1+shift),
3205     >            Horb)
3206         end do
3207
3208      else
3209
3210      if (.not.BA_push_get(mt_dcpl,npack1,'etmp',etmp(2),etmp(1)))
3211     > call errquit('psi_project_out_virtual1: out of stack',0,MA_ERR)
3212
3213         !call dcopy(2*npack1,0.0d0,0,dcpl_mb(etmp(1)),1)
3214         call Parallel_shared_vector_zero(.true.,
3215     >                                    2*npack1,dcpl_mb(etmp(1)))
3216         do i=1,neq(ms)
3217           call Pack_cc_dot(1,
3218     >            dcpl_mb(psi1(1) +(i-1)*npack1+shift),
3219     >            Horb,
3220     >            sum)
3221           call Pack_cc_daxpy(1,(-sum),
3222     >            dcpl_mb(psi1(1) +(i-1)*npack1+shift),
3223     >            dcpl_mb(etmp(1)))
3224         end do
3225         call D1dB_Vector_SumAll(2*npack1,dcpl_mb(etmp(1)))
3226         call daxpy_omp(2*npack1,1.0d0,dcpl_mb(etmp(1)),1,Horb,1)
3227
3228      if (.not.BA_pop_stack(etmp(2)))
3229     > call errquit('psi_project_out_virtual1:popping stack',0,MA_ERR)
3230
3231      end if
3232
3233      !**** project out virtual orbitals ****
3234      do jj=1,(kk-1)
3235        call Pack_cc_dot(1,
3236     >            dcpl_mb(psi1_excited(1) +(jj-1)*npack1+shifte),
3237     >            Horb,
3238     >            sum)
3239
3240        call Pack_cc_daxpy(1,(-sum),
3241     >            dcpl_mb(psi1_excited(1) +(jj-1)*npack1+shifte),
3242     >            Horb)
3243      end do
3244
3245
3246      return
3247      end
3248
3249
3250
3251
3252*     *******************************************
3253*     *				                *
3254*     *	         psi_project_out_virtual        *
3255*     *				                *
3256*     *******************************************
3257*
3258*    This routine projects out non-orthogonal components of Horb.
3259* This routine is needed for a KS minimizer.
3260*
3261      subroutine psi_project_out_virtual(ii,Horb)
3262      integer ii
3263      complex*16 Horb(*)
3264
3265#include "bafdecls.fh"
3266#include "errquit.fh"
3267#include "psi.fh"
3268
3269      integer ms,i,jj,kk,shift,shifte
3270      integer etmp(2)
3271      real*8  sum
3272
3273*     **** spin up orbital ****
3274      if (ii.le.ne_excited(1)) then
3275
3276         shift  = 0
3277         shifte = 0
3278         ms     = 1
3279         kk     = ii
3280*     **** spin down orbital ****
3281      else
3282         shift  = neq(1)*npack1
3283         shifte = ne_excited(1)*npack1
3284         ms     = 2
3285         kk     = ii-ne_excited(1)
3286      end if
3287
3288      !**** project out filled orbitals ****
3289      if (neq(ms).eq.ne(ms)) then
3290
3291         do i=1,ne(ms)
3292           call Pack_cc_dot(1,
3293     >            dcpl_mb(psi1(1) +(i-1)*npack1+shift),
3294     >            Horb,
3295     >            sum)
3296
3297           call Pack_cc_daxpy(1,(-sum),
3298     >            dcpl_mb(psi1(1) +(i-1)*npack1+shift),
3299     >            Horb)
3300         end do
3301
3302      else
3303      if (.not.BA_push_get(mt_dcpl,npack1,'etmp',etmp(2),etmp(1)))
3304     > call errquit('psi_project_out_virtual1: out of stack',0,MA_ERR)
3305
3306         !call dcopy(2*npack1,0.0d0,0,dcpl_mb(etmp(1)),1)
3307         call Parallel_shared_vector_zero(.true.,
3308     >                       2*npack1,dcpl_mb(etmp(1)))
3309         do i=1,neq(ms)
3310           call Pack_cc_dot(1,
3311     >            dcpl_mb(psi1(1) +(i-1)*npack1+shift),
3312     >            Horb,
3313     >            sum)
3314           call Pack_cc_daxpy(1,(-sum),
3315     >            dcpl_mb(psi1(1) +(i-1)*npack1+shift),
3316     >            dcpl_mb(etmp(1)))
3317         end do
3318         call D1dB_Vector_SumAll(2*npack1,dcpl_mb(etmp(1)))
3319         call daxpy(2*npack1,1.0d0,dcpl_mb(etmp(1)),1,Horb,1)
3320
3321      if (.not.BA_pop_stack(etmp(2)))
3322     > call errquit('psi_project_out_virtual:popping stack',0,MA_ERR)
3323
3324      end if
3325
3326      !**** project out virtual orbitals ****
3327      do jj=1,(kk)
3328        call Pack_cc_dot(1,
3329     >            dcpl_mb(psi1_excited(1) +(jj-1)*npack1+shifte),
3330     >            Horb,
3331     >            sum)
3332
3333        call Pack_cc_daxpy(1,(-sum),
3334     >            dcpl_mb(psi1_excited(1) +(jj-1)*npack1+shifte),
3335     >            Horb)
3336      end do
3337
3338
3339      return
3340      end
3341
3342
3343
3344************************ KS orbital Part ************************
3345
3346*
3347      subroutine psi_KS_update00(psi_number,
3348     >                         maxit_orbs,maxit_orb,
3349     >                         ks_algorithm,
3350     >                         precondition,
3351     >                         maxerror)
3352      implicit none
3353      integer psi_number
3354      integer maxit_orbs,maxit_orb
3355      integer ks_algorithm
3356      logical precondition
3357      real*8  maxerror
3358
3359#include "bafdecls.fh"
3360#include "psi.fh"
3361
3362*     **** local variables ****
3363      logical done
3364      integer i,j,neall
3365      real*8 error,error_out,tim1,tim2,tim,sum
3366
3367*     **** external functions ****
3368
3369      tim = 0.0d0
3370      neall = neq(1)+neq(2)
3371      j = 0
3372 2    j = j+1
3373        error = 0.0d0
3374        !do i=neall,1,-1
3375        do i=1,neall
3376
3377         !*** orthogonalize to lower orbitals  ****
3378         call psi_project_out_f_orb1(
3379     >           i,
3380     >           dcpl_mb(psi1(1)+(i-1)*npack1))
3381
3382         !*** normalize ****
3383         call Pack_cc_dot(1,
3384     >            dcpl_mb(psi1(1) +(i-1)*npack1),
3385     >            dcpl_mb(psi1(1) +(i-1)*npack1),
3386     >            sum)
3387         sum = 1.0d0/dsqrt(sum)
3388c         call Pack_c_SMul(1,sum,
3389c     >            dcpl_mb(psi1(1) +(i-1)*npack1),
3390c     >            dcpl_mb(psi1(1) +(i-1)*npack1))
3391         call Pack_c_SMul1(1,sum,
3392     >            dcpl_mb(psi1(1) +(i-1)*npack1))
3393
3394
3395
3396          if (ks_algorithm.eq.1) then
3397          call psi_KS_update_orb2(psi_number,precondition,maxit_orb,
3398     >                         maxerror,
3399     >                         0.1d0,i,error_out)
3400          else
3401          call psi_KS_update_orb(psi_number,precondition,maxit_orb,
3402     >                         maxerror,
3403     >                         0.1d0,i,error_out)
3404          end if
3405
3406          error = error+error_out
3407        end do
3408        error = error/dble(neall)
3409
3410        done = ((j.gt.maxit_orbs).or.(error.lt.maxerror))
3411      if (.not.done) go to 2
3412
3413      return
3414      end
3415
3416
3417
3418
3419*     ***********************************
3420*     *				        *
3421*     *	     psi_KS_update	        *
3422*     *				        *
3423*     ***********************************
3424
3425*    This routine (approximately) diagonalizes the KS matrix.
3426*
3427      subroutine psi_KS_update(psi_number,
3428     >                         ks_algorithm,
3429     >                         precondition,
3430     >                         maxerror)
3431      implicit none
3432      integer psi_number
3433      integer ks_algorithm
3434      logical precondition
3435      real*8 maxerror
3436
3437#include "bafdecls.fh"
3438#include "psi.fh"
3439
3440*     **** local variables ****
3441      logical done
3442      integer i,j,neall,maxit_orb,maxit_orbs
3443      real*8 error,error_out,tim1,tim2,tim,sum
3444
3445*     **** external functions ****
3446      integer  control_ks_maxit_orb,control_ks_maxit_orbs
3447      external control_ks_maxit_orb,control_ks_maxit_orbs
3448
3449      tim = 0.0d0
3450      neall = neq(1)+neq(2)
3451      maxit_orb  = control_ks_maxit_orb()   !*** should be read from rtdb ***
3452      maxit_orbs = control_ks_maxit_orbs()  !*** should be read from rtdb ***
3453      j = 0
3454 2    j = j+1
3455        error = 0.0d0
3456        !do i=neall,1,-1
3457        do i=1,neall
3458
3459         !*** orthogonalize to lower orbitals  ****
3460         call psi_project_out_f_orb1(
3461     >           i,
3462     >           dcpl_mb(psi1(1)+(i-1)*npack1))
3463
3464         !*** normalize ****
3465         call Pack_cc_dot(1,
3466     >            dcpl_mb(psi1(1) +(i-1)*npack1),
3467     >            dcpl_mb(psi1(1) +(i-1)*npack1),
3468     >            sum)
3469         sum = 1.0d0/dsqrt(sum)
3470c         call Pack_c_SMul(1,sum,
3471c     >            dcpl_mb(psi1(1) +(i-1)*npack1),
3472c     >            dcpl_mb(psi1(1) +(i-1)*npack1))
3473         call Pack_c_SMul1(1,sum,
3474     >            dcpl_mb(psi1(1) +(i-1)*npack1))
3475
3476
3477
3478          if (ks_algorithm.eq.1) then
3479          call psi_KS_update_orb2(psi_number,precondition,maxit_orb,
3480     >                         maxerror,
3481     >                         0.1d0,i,error_out)
3482          else
3483          call psi_KS_update_orb(psi_number,precondition,maxit_orb,
3484     >                         maxerror,
3485     >                         0.1d0,i,error_out)
3486          end if
3487
3488          error = error+error_out
3489        end do
3490        error = error/dble(neall)
3491
3492        done = ((j.gt.maxit_orbs).or.(error.lt.maxerror))
3493      if (.not.done) go to 2
3494
3495      return
3496      end
3497
3498
3499
3500
3501*     ***********************************
3502*     *                                 *
3503*     *      psi_KS_block_update        *
3504*     *                                 *
3505*     ***********************************
3506
3507*    This routine (approximately) diagonalizes the KS matrix.
3508*
3509*    input variables
3510*    - current_iteration
3511*    - maxerror
3512*    global shared variables output
3513*    - Enew
3514*    - deltae
3515*
3516      subroutine psi_KS_block_update(Enew,deltae,current_iteration,
3517     >                               maxit_orbs,maxerror)
3518      implicit none
3519      real*8  Enew,deltae
3520      integer current_iteration
3521      integer maxit_orbs
3522      real*8  maxerror
3523
3524#include "bafdecls.fh"
3525#include "errquit.fh"
3526#include "psi.fh"
3527
3528*     **** local variables ****
3529      real*8  deltat_min
3530      parameter (deltat_min=1.0d-2)
3531
3532      integer G0(2),S0(2)
3533      real*8  E0,dE0,Enew0
3534
3535      logical precondition
3536      common / cgsd_block2 / precondition
3537
3538      real*8 tmin,dte,sum0,sum1
3539      common / bfgs_block / tmin,dte,sum0,sum1
3540
3541      integer it,it_in
3542      real*8 deltat,tmin0,deltae0,deltac0
3543      real*8 max_sigma,dt,kappa
3544
3545      logical value
3546      integer neall
3547
3548*     **** external functions ****
3549      logical  control_precondition,psp_pawexist
3550      external control_precondition,psp_pawexist
3551      integer  control_lmbfgs_size
3552      external control_lmbfgs_size
3553      !integer  control_ks_maxit_orbs,psi_neq
3554      !external control_ks_maxit_orbs,psi_neq
3555      integer  psi_neq
3556      external psi_neq
3557      real*8   control_time_step
3558      external control_time_step
3559      real*8   psi_geodesic_energy0,psi_geodesic_denergy
3560      external psi_geodesic_energy0,psi_geodesic_denergy
3561      real*8   linesearch
3562      external linesearch
3563
3564      !maxit_orb  = control_ks_maxit_orb()   !*** should be read from rtdb ***
3565      !maxit_orbs = control_ks_maxit_orbs()  !*** should be read from rtdb ***
3566      it_in = maxit_orbs
3567
3568      dt = control_time_step()
3569
3570
3571!$OMP MASTER
3572      if (current_iteration.eq.1) then
3573         precondition = control_precondition()
3574         tmin  = 10*deltat_min
3575      end if
3576!$OMP END MASTER
3577
3578      neall = neq(1)+neq(2)
3579
3580*     **** allocate G0, S0 ****
3581      value = BA_push_get(mt_dcpl,npack1*neall,
3582     >                     'S0',S0(2),S0(1))
3583      value = value.and.
3584     >        BA_push_get(mt_dcpl,npack1*neall,
3585     >                     'G0',G0(2),G0(1))
3586      if (.not.value)
3587     >   call errquit('psi_KS_block_update:out of heap memory',0,MA_ERR)
3588
3589      call Parallel_shared_vector_zero(.true.,2*npack1*neall,
3590     >                                 dcpl_mb(G0(1)))
3591
3592*     ***** get the initial direction ****
3593      if (pawexist) then
3594         call psi_1get_STgradient(dcpl_mb(S0(1)),dcpl_mb(G0(1)),Enew0)
3595      else
3596         if (precondition) then
3597            call psi_1get_TMgradient0(dcpl_mb(G0(1)),Enew0)
3598         else
3599            call psi_1get_Tgradient0(dcpl_mb(G0(1)),Enew0)
3600         end if
3601      end if
3602      E0 = Enew0
3603
3604*     ***** use the initial gradient for the direction ****
3605      call pspw_lmbfgs_init(control_lmbfgs_size(),dcpl_mb(G0(1)))
3606      call Grsm_gg_Copy(npack1,neall,
3607     >                  dcpl_mb(G0(1)),
3608     >                  dcpl_mb(S0(1)))
3609
3610      do it=2,it_in
3611*        **** initialize the geoedesic line data structure ****
3612         call geodesic_start(dcpl_mb(S0(1)),max_sigma,dE0)
3613
3614*        ******* line search *********
3615         if ((tmin.gt.deltat_min).and.(tmin.lt.1.0d4)) then
3616            deltat = tmin
3617         else
3618            deltat = deltat_min
3619         end if
3620 20      continue
3621c         write(*,*) "     --- tmin0=",tmin
3622         tmin0 = tmin
3623         deltae0 = deltae
3624         Enew0 = linesearch(0.0d0,E0,dE0,deltat,
3625     >                        psi_geodesic_energy0,
3626     >                        psi_geodesic_denergy,
3627     >                        0.50d0,tmin0,deltae0,2)
3628
3629!$OMP MASTER
3630         tmin = tmin0
3631         deltae = deltae0
3632         Enew = Enew0
3633!$OMP END MASTER
3634!$OMP BARRIER
3635c         write(*,*) "     --- tmin=",tmin
3636c         write(*,*) "     --- it,Enew,deltae=",it,Enew,deltae
3637
3638         call psi_geodesic_final(tmin)
3639
3640*        **** exit loop early ****
3641         if (dabs(deltae).lt.maxerror) then
3642            if (.not.precondition) go to 30
3643            precondition = .false.
3644         end if
3645
3646*        **** get the new gradient ****
3647         if (pawexist) then
3648              call psi_2get_STgradient(2,dcpl_mb(S0(1)),
3649     >                                   dcpl_mb(G0(1)),Enew0)
3650         else
3651            if (precondition) then
3652               call psi_2get_TMgradient0(2,dcpl_mb(G0(1)),Enew0)
3653            else
3654               call psi_2get_Tgradient0(2,dcpl_mb(G0(1)),Enew0)
3655            end if
3656         end if
3657!$OMP MASTER
3658      E0 = Enew0
3659!$OMP END MASTER
3660!$OMP BARRIER
3661
3662         call pspw_lmbfgs(tmin,dcpl_mb(G0(1)),dcpl_mb(S0(1)))
3663         call psi_2to1()
3664
3665      end do
3666
3667*     **** initialize the geoedesic line data structure ****
3668      call geodesic_start(dcpl_mb(S0(1)),max_sigma,dE0)
3669
3670*     ******* line search *********
3671      if ((tmin.gt.deltat_min).and.(tmin.lt.1.0d4)) then
3672         deltat = tmin
3673      else
3674         deltat = deltat_min
3675      end if
3676
3677 25   continue
3678      tmin0 = tmin
3679      deltae0 = deltae
3680      Enew0 = linesearch(0.0d0,E0,dE0,deltat,
3681     >                        psi_geodesic_energy0,
3682     >                        psi_geodesic_denergy,
3683     >                        0.50d0,tmin0,deltae0,2)
3684!$OMP MASTER
3685         tmin = tmin0
3686         deltae = deltae0
3687         Enew = Enew0
3688!$OMP END MASTER
3689!$OMP BARRIER
3690
3691      call psi_geodesic_final(tmin)
3692
3693 30   CONTINUE
3694
3695
3696*     **** free memory ****
3697      value = BA_pop_stack(G0(2))
3698      value = value.and.
3699     >        BA_pop_stack(S0(2))
3700      if (.not. value)
3701     >  call errquit('bfgsminimize:error freeing heap memory',0, MA_ERR)
3702
3703      call psi_2to1()
3704
3705      return
3706      end
3707
3708
3709
3710*     ***********************************
3711*     *				        *
3712*     *	     psi_KS_update_orb	        *
3713*     *				        *
3714*     ***********************************
3715
3716*    This routine performs a KS update on orbital i
3717*
3718      subroutine psi_KS_update_orb(psi_number,
3719     >                             precondition,maxiteration,
3720     >                             maxerror,perror,i,
3721     >                             error_out)
3722      implicit none
3723      integer psi_number
3724      logical precondition
3725      integer maxiteration
3726      real*8  maxerror,perror
3727      integer i
3728      real*8 error_out
3729
3730#include "bafdecls.fh"
3731#include "psi.fh"
3732#include "errquit.fh"
3733
3734*     **** local variables ****
3735      integer MASTER,taskid
3736      parameter (MASTER=0)
3737
3738      logical value,done,oneloop
3739      integer it
3740      real*8 e0,eold,error0,de0,lmbda_r0,lmbda_r1
3741      real*8 theta,sigma
3742      integer r1(2),t0(2),t(2),g(2)
3743      integer psi_ptr
3744
3745      if (psi_number.eq.1) then
3746         psi_ptr=psi1(1)
3747      else
3748         psi_ptr=psi2(1)
3749      end if
3750
3751      call Parallel_taskid(taskid)
3752
3753      value = BA_push_get(mt_dcpl,npack1,'t0',t0(2),t0(1))
3754      value = value.and.
3755     >        BA_push_get(mt_dcpl,npack1,'r1',r1(2),r1(1))
3756      value = value.and.
3757     >        BA_push_get(mt_dcpl,npack1,'g',g(2),g(1))
3758      value = value.and.
3759     >        BA_push_get(mt_dcpl,npack1,'t',t(2),t(1))
3760      if (.not. value) call errquit(
3761     >     'psi_KS_update_orb: out of stack memory',0, MA_ERR)
3762
3763      done = .false.
3764      error0 = 0.0d0
3765      e0 = 0.0d0
3766      theta = -3.14159d0/600.0d0
3767      lmbda_r0 = 1.0d0
3768      it = 0
3769 2    continue
3770
3771         it = it + 1
3772         eold = e0
3773
3774*        *** calculate residual (steepest descent) direction for a single band ***
3775         call psi_get_gradient_orb(psi_number,i,dcpl_mb(g(1)))
3776         call Pack_cc_dot(1,dcpl_mb(psi_ptr+(i-1)*npack1),
3777     >                   dcpl_mb(g(1)),
3778     >                    e0)
3779         !dbl_mb(eig(1)+i-1) = e0
3780         e0 = -e0
3781
3782         done = ((it.gt.maxiteration)
3783     >           .or.
3784     >           (dabs(e0-eold).lt.maxerror))
3785
3786         if (done) go to 4
3787
3788*        **** preconditioning ****
3789         if (precondition) then
3790           call ke_Precondition(npack1,1,
3791     >                     dcpl_mb(psi_ptr+(i-1)*npack1),
3792     >                     dcpl_mb(g(1)))
3793
3794         end if
3795
3796c        call Pack_c_Copy(1,dcpl_mb(g(1)),dcpl_mb(r1(1)))
3797c        call Pack_cc_daxpy(1,(e0),
3798c    >                 dcpl_mb(psi_ptr+(i-1)*npack1),
3799c    >                 dcpl_mb(r1(1)))
3800         call Pack_c_Copy(1,dcpl_mb(g(1)),dcpl_mb(r1(1)))
3801         call psi_project_out_orb(psi_number,i,dcpl_mb(r1(1)))
3802
3803
3804
3805
3806*        *** determine conjuagate direction ***
3807         call Pack_cc_dot(1,dcpl_mb(r1(1)),
3808     >                   dcpl_mb(r1(1)),
3809     >                   lmbda_r1)
3810         call Pack_c_Copy(1,dcpl_mb(r1(1)),dcpl_mb(t(1)))
3811
3812         if (it.gt.1) then
3813         call Pack_cc_daxpy(1,(lmbda_r1/lmbda_r0),
3814     >                   dcpl_mb(t0(1)),
3815     >                   dcpl_mb(t(1)))
3816         end if
3817         lmbda_r0 = lmbda_r1
3818         oneloop = .true.
3819 3       call Pack_c_Copy(1,dcpl_mb(t(1)),dcpl_mb(t0(1)))
3820
3821
3822
3823*        **** project out psi components from t - may not be needed! ****
3824         !call psi_project_out_orb(psi_number,i,dcpl_mb(t(1)))
3825c!        call psi_project_out_orb(psi_number,i,dcpl_mb(t(1)))
3826c!        call Pack_cc_dot(1,dcpl_mb(psi_ptr+(i-1)*npack1),
3827c!    >                   dcpl_mb(t(1)),
3828c!    >                    de0)
3829c!        de0 = -de0
3830c!        call Pack_cc_daxpy(1,(de0),
3831c!    >                 dcpl_mb(psi_ptr+(i-1)*npack1),
3832c!    >                 dcpl_mb(t(1)))
3833
3834
3835*        *** normalize search direction, t ****
3836         call Pack_cc_dot(1,dcpl_mb(t(1)),
3837     >                   dcpl_mb(t(1)),
3838     >                   sigma)
3839         sigma = dsqrt(sigma)
3840         de0 = 1.0d0/sigma
3841c         call Pack_c_SMul(1,de0,dcpl_mb(t(1)),dcpl_mb(t(1)))
3842         call Pack_c_SMul1(1,de0,dcpl_mb(t(1)))
3843
3844
3845*        **** compute de0 = <t|g> ****
3846         call Pack_cc_dot(1,dcpl_mb(t(1)),
3847     >                   dcpl_mb(g(1)),
3848     >                   de0)
3849
3850*        *** bad direction ***
3851         if ((de0.lt.0.0d0).and.oneloop) then
3852           call Pack_c_Copy(1,dcpl_mb(g(1)),dcpl_mb(t(1)))
3853           oneloop = .false.
3854           go to 3
3855         end if
3856
3857         de0 = -2.0d0*de0
3858         call psi_linesearch_update2(psi_number,i,
3859     >                              theta,e0,de0,
3860     >                              dcpl_mb(t(1)),
3861     >                              sigma,
3862     >                              dcpl_mb(t0(1)))
3863
3864      go to 2
3865
3866
3867*     **** release stack memory ****
3868 4    value =           BA_pop_stack(t(2))
3869      value = value.and.BA_pop_stack(g(2))
3870      value = value.and.BA_pop_stack(r1(2))
3871      value = value.and.BA_pop_stack(t0(2))
3872      if (.not. value) call errquit(
3873     >     'psi_KS_update_orb: popping stack memory',1, MA_ERR)
3874
3875c      write(*,*) "iterations=",it," eig=",e0," error=",error_out,
3876c     >           theta
3877      error_out = dabs(e0-eold)
3878      return
3879      end
3880
3881
3882
3883
3884*     ***********************************
3885*     *				        *
3886*     *	     psi_KS_update_orb2	        *
3887*     *				        *
3888*     ***********************************
3889
3890*    This routine performs a RMM-DIIS KS update on orbital i
3891*
3892      subroutine psi_KS_update_orb2(psi_number,
3893     >                             precondition,maxiteration,
3894     >                             maxerror,perror,i,
3895     >                             error_out)
3896      implicit none
3897      integer psi_number
3898      logical precondition
3899      integer maxiteration
3900      real*8  maxerror,perror
3901      integer i
3902      real*8 error_out
3903
3904#include "bafdecls.fh"
3905#include "errquit.fh"
3906#include "psi.fh"
3907
3908*     **** local variables ****
3909      integer MASTER,taskid
3910      parameter (MASTER=0)
3911
3912      logical value,done
3913      integer it
3914      real*8 sigma,e0,eold,error0
3915      real*8 lambda
3916      integer r1(2),g1(2)
3917      integer psi_ptr
3918
3919      if (psi_number.eq.1) then
3920         psi_ptr=psi1(1)
3921      else
3922         psi_ptr=psi2(1)
3923      end if
3924
3925      call Parallel_taskid(taskid)
3926
3927*     **** allocate memory ****
3928      value =            BA_push_get(mt_dcpl,npack1,'r1',r1(2),r1(1))
3929      value = value.and. BA_push_get(mt_dcpl,npack1,'g1',g1(2),g1(1))
3930      if (.not. value) call errquit(
3931     >     'psi_KS_update_orb2: out of stack memory',0, MA_ERR)
3932
3933*     **** set lambda ***
3934      lambda = 0.1d0
3935
3936
3937*     *** calculate residual (steepest descent) direction for a single band ***
3938      call psi_get_gradient_orb(psi_number,i,dcpl_mb(g1(1)))
3939      call Pack_cc_dot(1,
3940     >                 dcpl_mb(psi_ptr+(i-1)*npack1),
3941     >                 dcpl_mb(g1(1)),
3942     >                 e0)
3943      call Pack_cc_dot(1,
3944     >                 dcpl_mb(psi_ptr+(i-1)*npack1),
3945     >                 dcpl_mb(psi_ptr+(i-1)*npack1),
3946     >                 sigma)
3947      e0 = e0/sigma
3948      call Pack_c_SMul(1,e0,
3949     >                 dcpl_mb(psi_ptr+(i-1)*npack1),
3950     >                 dcpl_mb(r1(1)))
3951      call Pack_cc_daxpy(1,(-1.0d0),
3952     >                 dcpl_mb(g1(1)),
3953     >                 dcpl_mb(r1(1)))
3954
3955        !write(*,*) "i=",i,"it=",0, " eig=",e0,sigma
3956
3957*     ***** rmmdiis start *****
3958      call pspw_rmmdiis_start(lambda,
3959     >                        dcpl_mb(r1(1)),
3960     >                        dcpl_mb(psi_ptr+(i-1)*npack1))
3961
3962
3963      done = .false.
3964      error0 = 0.0d0
3965      it = 0
3966 2    continue
3967
3968         it = it + 1
3969         eold = e0
3970
3971*        *** calculate residual (steepest descent) direction for a single band ***
3972         call psi_get_gradient_orb(psi_number,i,dcpl_mb(g1(1)))
3973         call Pack_cc_dot(1,
3974     >                    dcpl_mb(psi_ptr+(i-1)*npack1),
3975     >                    dcpl_mb(g1(1)),
3976     >                    e0)
3977         call Pack_cc_dot(1,
3978     >                    dcpl_mb(psi_ptr+(i-1)*npack1),
3979     >                    dcpl_mb(psi_ptr+(i-1)*npack1),
3980     >                    sigma)
3981         e0 = e0/sigma
3982         call Pack_c_SMul(1,(e0),
3983     >                    dcpl_mb(psi_ptr+(i-1)*npack1),
3984     >                    dcpl_mb(r1(1)))
3985         call Pack_cc_daxpy(1,(-1.0d0),
3986     >                    dcpl_mb(g1(1)),
3987     >                    dcpl_mb(r1(1)))
3988         !e0 = -e0
3989        !write(*,*) "i=",i,"it=",it, " eig=",e0,sigma,dabs(e0-eold)
3990
3991         done = ((it.gt.maxiteration)
3992     >           .or.
3993     >           (dabs(e0-eold).lt.maxerror))
3994
3995         if (done) go to 4
3996
3997*        ***** rmmdiis update *****
3998         call pspw_rmmdiis(lambda,
3999     >                     dcpl_mb(r1(1)),
4000     >                     dcpl_mb(psi_ptr+(i-1)*npack1))
4001
4002      go to 2
4003
4004*     ***** normalize psi ****
4005 4    call Pack_cc_dot(1,
4006     >                 dcpl_mb(psi_ptr+(i-1)*npack1),
4007     >                 dcpl_mb(psi_ptr+(i-1)*npack1),
4008     >                 sigma)
4009c      call Pack_c_SMul(1,(1.0d0/dsqrt(sigma)),
4010c     >                 dcpl_mb(psi_ptr+(i-1)*npack1),
4011c     >                 dcpl_mb(psi_ptr+(i-1)*npack1))
4012      call Pack_c_SMul1(1,(1.0d0/dsqrt(sigma)),
4013     >                 dcpl_mb(psi_ptr+(i-1)*npack1))
4014
4015
4016*     **** release stack memory ****
4017      value =           BA_pop_stack(g1(2))
4018      value = value.and.BA_pop_stack(r1(2))
4019      if (.not. value) call errquit(
4020     >     'psi_KS_update_orb2: popping stack memory',1, MA_ERR)
4021      error_out = dabs(e0-eold)
4022
4023c       write(*,*) "i=",i,"iterations=",it," eig=",e0,
4024c     >            " error=",error_out,
4025c     >            lambda
4026      return
4027      end
4028
4029
4030
4031
4032
4033
4034
4035
4036*     ***********************************
4037*     *				        *
4038*     *	     psi_linesearch_update	*
4039*     *				        *
4040*     ***********************************
4041
4042*    This routine performs a linesearch on orbital i, in the direction t.
4043* This routine is needed for a KS minimizer.
4044*  e0 = <orb|g>
4045*  de0 = 2*<t|g>
4046*
4047      subroutine psi_linesearch_update(psi_number,i,theta,e0,de0,t)
4048      implicit none
4049#include "errquit.fh"
4050      integer psi_number
4051      integer i
4052      real*8  theta
4053      real*8  e0,de0
4054      complex*16 t(*) !search direction
4055
4056#include "bafdecls.fh"
4057#include "psi.fh"
4058
4059*     **** local variables ****
4060      logical value
4061      integer orb(2),g(2),psi_ptr
4062      real*8 x,y,pi,e1
4063
4064      if (psi_number.eq.1) then
4065         psi_ptr=psi1(1)
4066      else
4067         psi_ptr=psi2(1)
4068      end if
4069
4070      pi = 4.0d0*datan(1.0d0)
4071
4072*     **** allocate stack memory ****
4073      value = BA_push_get(mt_dcpl,npack1,'orb',
4074     >                       orb(2),orb(1))
4075      value = value.and.
4076     >        BA_push_get(mt_dcpl,npack1,'g',
4077     >                       g(2),g(1))
4078      if (.not. value) call errquit(
4079     >     'psi_linesearch_update: out of stack memory',0, MA_ERR)
4080
4081
4082      call Pack_c_Copy(1,dcpl_mb(psi_ptr+(i-1)*npack1),
4083     >                   dcpl_mb(orb(1)))
4084
4085*     **** orb2 = orb*cos(pi/300) + t*sin(pi/300) ****
4086      !theta = pi/300.0d0
4087      x = cos(theta)
4088      y = sin(theta)
4089      call Pack_c_SMul(1,x,
4090     >                  dcpl_mb(orb(1)),
4091     >                  dcpl_mb(psi_ptr+(i-1)*npack1))
4092      call Pack_cc_daxpy(1,y,
4093     >                   t,
4094     >                   dcpl_mb(psi_ptr+(i-1)*npack1))
4095
4096*     *** determine theta ***
4097      call psi_get_gradient_orb(psi_number,i,dcpl_mb(g(1)))
4098
4099      call Pack_cc_dot(1,dcpl_mb(psi_ptr+(i-1)*npack1),
4100     >                   dcpl_mb(g(1)),
4101     >                   e1)
4102      e1 = -e1
4103      x = (e0 - e1 + 0.5d0*de0*sin(2*theta))
4104     >    /(1.0d0-cos(2*theta))
4105      theta = 0.5d0*datan(0.5d0*de0/x)
4106
4107c     call Pack_cc_dot(1,t,
4108c    >                 dcpl_mb(g(1)),
4109c    >                 de1)
4110c     de1 = -2.0d0*de1
4111c     theta  = -de1*(pi/300.0d0)/(de1-de0)
4112
4113      !write(*,*) "i,theta,e1:",i,theta,e1
4114
4115
4116*     **** orb2 = orb*cos(theta) + t*sin(theta) ****
4117      x = cos(theta)
4118      y = sin(theta)
4119      call Pack_c_SMul(1,x,
4120     >                  dcpl_mb(orb(1)),
4121     >                  dcpl_mb(psi_ptr+(i-1)*npack1))
4122      call Pack_cc_daxpy(1,y,
4123     >                   t,
4124     >                   dcpl_mb(psi_ptr+(i-1)*npack1))
4125
4126*     **** update orb2_r and H*orb2 ****
4127      !call electron_run_orb(i,dcpl_mb(psi_ptr))
4128c     call psi_get_gradient_orb(psi_number,i,dcpl_mb(g(1)))
4129c     call Pack_cc_dot(1,dcpl_mb(psi_ptr+(i-1)*npack1),
4130c    >                   dcpl_mb(g(1)),
4131c    >                   e2)
4132c     e2 = -e2
4133c     call Pack_cc_dot(1,t,
4134c    >                 dcpl_mb(g(1)),
4135c    >                 de2)
4136c     de2 = -2.0d0*de2
4137
4138c     write(*,*) "i,theta,es:",i,theta,e0,e1,e2
4139c     write(*,*)
4140
4141*     **** release stack memory ****
4142      value =           BA_pop_stack(g(2))
4143      value = value.and.BA_pop_stack(orb(2))
4144      if (.not. value) call errquit(
4145     >     'psi_linesearch_update: popping stack memory',1, MA_ERR)
4146
4147      return
4148      end
4149
4150*     ***********************************
4151*     *				        *
4152*     *	     psi_linesearch_update2	*
4153*     *				        *
4154*     ***********************************
4155
4156*    This routine performs a linesearch on orbital i, in the direction t.
4157* This routine is needed for a KS minimizer.
4158*  e0 = <orb|g>
4159*  de0 = 2*<t|g>
4160*
4161      subroutine psi_linesearch_update2(psi_number,i,theta,e0,de0,t,
4162     >                                  sigma,tau_t)
4163      implicit none
4164#include "errquit.fh"
4165      integer psi_number
4166      integer i
4167      real*8  theta
4168      real*8  e0,de0
4169      complex*16 t(*)     !search direction
4170
4171      real*8     sigma
4172      complex*16 tau_t(*) !parallel transported search direction
4173
4174#include "bafdecls.fh"
4175#include "psi.fh"
4176
4177*     **** local variables ****
4178      logical value
4179      integer orb(2),g(2),psi_ptr
4180      real*8 x,y,pi,e1
4181
4182      if (psi_number.eq.1) then
4183         psi_ptr=psi1(1)
4184      else
4185         psi_ptr=psi2(1)
4186      end if
4187
4188      pi = 4.0d0*datan(1.0d0)
4189
4190*     **** allocate stack memory ****
4191      value = BA_push_get(mt_dcpl,npack1,'orb',
4192     >                       orb(2),orb(1))
4193      value = value.and.
4194     >        BA_push_get(mt_dcpl,npack1,'g',
4195     >                       g(2),g(1))
4196      if (.not. value) call errquit(
4197     >     'psi_linesearch_update: out of stack memory',0, MA_ERR)
4198
4199
4200      call Pack_c_Copy(1,dcpl_mb(psi_ptr+(i-1)*npack1),
4201     >                   dcpl_mb(orb(1)))
4202
4203*     **** orb2 = orb*cos(pi/300) + t*sin(pi/300) ****
4204      !theta = pi/300.0d0
4205      x = cos(theta)
4206      y = sin(theta)
4207      call Pack_c_SMul(1,x,
4208     >                  dcpl_mb(orb(1)),
4209     >                  dcpl_mb(psi_ptr+(i-1)*npack1))
4210      call Pack_cc_daxpy(1,y,
4211     >                   t,
4212     >                   dcpl_mb(psi_ptr+(i-1)*npack1))
4213
4214*     *** determine theta ***
4215      call psi_get_gradient_orb(psi_number,i,dcpl_mb(g(1)))
4216
4217      call Pack_cc_dot(1,dcpl_mb(psi_ptr+(i-1)*npack1),
4218     >                   dcpl_mb(g(1)),
4219     >                   e1)
4220      e1 = -e1
4221      x = (e0 - e1 + 0.5d0*de0*sin(2*theta))
4222     >    /(1.0d0-cos(2*theta))
4223      theta = 0.5d0*datan(0.5d0*de0/x)
4224
4225      x = cos(theta)
4226      y = sin(theta)
4227
4228*     **** tau_t = (-orb*sin(theta) + t*cos(theta))*sigma ****
4229      call Pack_c_SMul(1,(-y),
4230     >                  dcpl_mb(orb(1)),
4231     >                  tau_t)
4232      call Pack_cc_daxpy(1,x,
4233     >                   t,
4234     >                   tau_t)
4235c      call Pack_c_SMul(1,sigma,
4236c     >                  tau_t,
4237c     >                  tau_t)
4238      call Pack_c_SMul1(1,sigma,tau_t)
4239
4240*     **** orb2 = orb*cos(theta) + t*sin(theta) ****
4241      call Pack_c_SMul(1,x,
4242     >                  dcpl_mb(orb(1)),
4243     >                  dcpl_mb(psi_ptr+(i-1)*npack1))
4244      call Pack_cc_daxpy(1,y,
4245     >                   t,
4246     >                   dcpl_mb(psi_ptr+(i-1)*npack1))
4247
4248
4249*     **** release stack memory ****
4250      value =           BA_pop_stack(g(2))
4251      value = value.and.BA_pop_stack(orb(2))
4252      if (.not. value) call errquit(
4253     >     'psi_linesearch_update: popping stack memory',1, MA_ERR)
4254
4255      return
4256      end
4257
4258
4259
4260*     ***************************
4261*     *				*
4262*     *	     psi_set_orb	*
4263*     *				*
4264*     ***************************
4265
4266*    This routine copies an orbital, orb, into the ith psi of psi1.
4267* This routine is needed for a KS minimizer.
4268*
4269      subroutine psi_set_orb(psi_number,i,orb)
4270      implicit none
4271      integer psi_number
4272      integer i
4273      complex*16 orb(*)
4274
4275#include "bafdecls.fh"
4276#include "psi.fh"
4277
4278*     **** local variables ****
4279      integer index,psi_ptr
4280
4281      if (psi_number.eq.1) then
4282         psi_ptr=psi1(1)
4283      else
4284         psi_ptr=psi2(1)
4285      end if
4286
4287      index = (i-1)*npack1
4288
4289      call zcopy(npack1,
4290     >           orb, 1,
4291     >           dcpl_mb(psi_ptr+index),1)
4292      return
4293      end
4294
4295
4296*     ***************************
4297*     *				*
4298*     *	     psi_get_orb	*
4299*     *				*
4300*     ***************************
4301
4302*    This routine copies the ith psi of psi1 into an orbital, orb.
4303* This routine is needed for a KS minimizer.
4304*
4305      subroutine psi_get_orb(psi_number,i,orb)
4306      implicit none
4307      integer psi_number
4308      integer i
4309      complex*16 orb(*)
4310
4311#include "bafdecls.fh"
4312#include "psi.fh"
4313
4314*     **** local variables ****
4315      integer index,psi_ptr
4316
4317
4318      if (psi_number.eq.1) then
4319         psi_ptr=psi1(1)
4320      else
4321         psi_ptr=psi2(1)
4322      end if
4323
4324      index = (i-1)*npack1
4325
4326      call zcopy(npack1,
4327     >           dcpl_mb(psi_ptr+index), 1,
4328     >           orb, 1)
4329      return
4330      end
4331
4332*     ***********************************
4333*     *				        *
4334*     *	     psi_get_gradient_orb	*
4335*     *				        *
4336*     ***********************************
4337
4338*    This routine returns the Hpsi(i).
4339* This routine is needed for a KS minimizer.
4340*
4341      subroutine psi_get_gradient_orb(psi_number,i,Horb)
4342      implicit none
4343      integer psi_number
4344      integer i
4345      complex*16 Horb(*)
4346
4347#include "bafdecls.fh"
4348#include "psi.fh"
4349
4350*     **** local variables ****
4351      integer psi_ptr
4352
4353      if (psi_number.eq.1) then
4354         psi_ptr=psi1(1)
4355      else
4356         psi_ptr=psi2(1)
4357      end if
4358
4359      call electron_run_orb(i,dcpl_mb(psi_ptr))
4360      call electron_get_gradient_orb(i,Horb)
4361
4362      return
4363      end
4364
4365
4366*     *******************************************
4367*     *				                *
4368*     *	         psi_project_out_orb           *
4369*     *				                *
4370*     *******************************************
4371*
4372*    This routine projects out non-orthogonal components of Horb.
4373* This routine is needed for a KS minimizer.
4374*
4375      subroutine psi_project_out_orb(psi_number,i,Horb)
4376      implicit none
4377#include "errquit.fh"
4378      integer psi_number
4379      integer i
4380      complex*16 Horb(*)
4381
4382#include "bafdecls.fh"
4383#include "psi.fh"
4384
4385*     **** local variables ****
4386      logical ok
4387      integer ii,n,psi_ptr,np
4388      integer x(2)
4389      real*8  sum
4390
4391      call Parallel_np(np)
4392
4393*     **** allocate stack memory ****
4394      ok = BA_push_get(mt_dbl,ne(1),'x',x(2),x(1))
4395      if (.not.ok)
4396     > call errquit('psi_project_out_orb: out of stack memory',0,
4397     &       MA_ERR)
4398
4399
4400      if (psi_number.eq.1) then
4401         psi_ptr=psi1(1)
4402      else
4403         psi_ptr=psi2(1)
4404      end if
4405
4406*     **** spin up orbital ****
4407      if (i.le.ne(1)) then
4408
4409        ii = i
4410!       do n=1,(ii)
4411!          call Pack_cc_dot(1,
4412!    >            dcpl_mb(psi_ptr +(n-1)*npack1),
4413!    >            Horb,
4414!    >            sum)
4415!          call daxpy(2*npack1,
4416!    >               (-sum),
4417!    >               dcpl_mb(psi_ptr+(n-1)*npack1),1,
4418!    >               Horb,1)
4419!       end do
4420        call Pack_cc_ndot(1,ii,
4421     >            dcpl_mb(psi_ptr),
4422     >            Horb,
4423     >            dbl_mb(x(1)))
4424        do n=1,(ii)
4425           call daxpy(2*npack1,
4426     >               (-dbl_mb(x(1)+n-1)),
4427     >               dcpl_mb(psi_ptr+(n-1)*npack1),1,
4428     >               Horb,1)
4429        end do
4430
4431
4432
4433*     **** spin down orbital ****
4434      else
4435
4436
4437        ii = i - ne(1)
4438        do n=(ne(1)+1),(ne(1)+ii)
4439           call Pack_cc_dot(1,
4440     >            dcpl_mb(psi_ptr +(n-1)*npack1),
4441     >            Horb,
4442     >            sum)
4443           call daxpy(2*npack1,
4444     >               (-sum),
4445     >               dcpl_mb(psi_ptr+(n-1)*npack1),1,
4446     >               Horb,1)
4447        end do
4448
4449
4450      end if
4451
4452*     **** release stack memory ****
4453      ok = BA_pop_stack(x(2))
4454      if (.not. ok)
4455     > call errquit('psi_project_out_orb: poping stack memory',0,
4456     &       MA_ERR)
4457
4458      return
4459      end
4460
4461
4462
4463
4464
4465
4466*     ***************************
4467*     *				*
4468*     *	     psi_set_density	*
4469*     *				*
4470*     ***************************
4471
4472*    This routine sets the densities and potentials in psi and electron.
4473* This routine is needed for a band by band minimizer.
4474*
4475      subroutine psi_set_density(psi_number,rho)
4476      implicit none
4477      integer psi_number
4478      real*8 rho(*)
4479
4480
4481#include "bafdecls.fh"
4482#include "psi.fh"
4483
4484*     ***** rhoall common block ****
4485      integer rho1_all(2)
4486      integer rho2_all(2)
4487      common / rhoall_block / rho1_all,rho2_all
4488
4489*     **** local variables ****
4490      integer rho_ptr,dng_ptr,rho_all_ptr
4491
4492      if (psi_number.eq.1) then
4493        rho_ptr     = rho1(1)
4494        dng_ptr     = dng1(1)
4495        rho_all_ptr = rho1_all(1)
4496      else
4497        rho_ptr     = rho2(1)
4498        dng_ptr     = dng2(1)
4499        rho_all_ptr = rho2_all(1)
4500      end if
4501
4502c      call dcopy(4*nfft3d,
4503c     >           rho, 1,
4504c     >           dbl_mb(rho_ptr),1)
4505
4506      call Parallel_shared_vector_copy(.true.,4*nfft3d,
4507     >                                 rho,dbl_mb(rho_ptr))
4508
4509      call electron_gen_dng_dnall(dbl_mb(rho_ptr),
4510     >                            dcpl_mb(dng_ptr),
4511     >                            dbl_mb(rho_all_ptr))
4512      call electron_gen_scf_potentials(dbl_mb(rho_ptr),
4513     >                            dcpl_mb(dng_ptr),
4514     >                            dbl_mb(rho_all_ptr))
4515      call electron_gen_vall()
4516      return
4517      end
4518
4519
4520*     ***************************
4521*     *				*
4522*     *	     psi_get_density	*
4523*     *				*
4524*     ***************************
4525
4526*    This routine gets the densities in psi.
4527* This routine is needed for a band by band minimizer.
4528*
4529      subroutine psi_get_density(psi_number,rho)
4530      implicit none
4531      integer psi_number
4532      real*8 rho(*)
4533
4534
4535#include "bafdecls.fh"
4536#include "psi.fh"
4537
4538*     ***** rhoall common block ****
4539      integer rho1_all(2)
4540      integer rho2_all(2)
4541      common / rhoall_block / rho1_all,rho2_all
4542
4543*     **** local variables ****
4544      integer rho_ptr
4545
4546      if (psi_number.eq.1) then
4547        rho_ptr = rho1(1)
4548      else
4549        rho_ptr = rho2(1)
4550      end if
4551
4552c      call dcopy(4*nfft3d,
4553c     >           dbl_mb(rho_ptr),1,
4554c     >           rho,1)
4555      call Parallel_shared_vector_copy(.true.,4*nfft3d,
4556     >                                 dbl_mb(rho_ptr),rho)
4557      return
4558      end
4559
4560*     ***************************
4561*     *                         *
4562*     *      psi_write_density  *
4563*     *                         *
4564*     ***************************
4565
4566*    This routine writes the densities in psi to disk.
4567* This routine is needed for a band by band minimizer.
4568*
4569      subroutine psi_write_density(psi_number)
4570      implicit none
4571      integer psi_number
4572
4573#include "bafdecls.fh"
4574#include "psi.fh"
4575
4576*     ***** rhoall common block ****
4577      integer rho1_all(2)
4578      integer rho2_all(2)
4579      common / rhoall_block / rho1_all,rho2_all
4580
4581*     **** local variables ****
4582      integer rho_ptr
4583
4584      if (psi_number.eq.1) then
4585        rho_ptr = rho1(1)
4586      else
4587        rho_ptr = rho2(1)
4588      end if
4589
4590      call rho_write(ispin,dbl_mb(rho_ptr))
4591      return
4592      end
4593
4594*     ***************************
4595*     *                         *
4596*     *   psi_try_read_density  *
4597*     *                         *
4598*     ***************************
4599
4600*    This routine reads the densities from disk to psi.
4601* This routine is needed for a band by band minimizer.
4602*
4603      logical function psi_try_read_density(psi_number)
4604      implicit none
4605      integer psi_number
4606
4607#include "bafdecls.fh"
4608#include "psi.fh"
4609
4610*     ***** rhoall common block ****
4611      integer rho1_all(2)
4612      integer rho2_all(2)
4613      common / rhoall_block / rho1_all,rho2_all
4614
4615*     **** local variables ****
4616      logical ok
4617      integer rho_ptr
4618
4619*     **** external functions ****
4620      logical  rho_check_header
4621      external rho_check_header
4622
4623      if (psi_number.eq.1) then
4624        rho_ptr = rho1(1)
4625      else
4626        rho_ptr = rho2(1)
4627      end if
4628
4629      if (rho_check_header(ispin,.false.)) then
4630         call rho_read(ispin,dbl_mb(rho_ptr))
4631         ok = .true.
4632      else
4633         ok = .false.
4634      end if
4635
4636      psi_try_read_density = ok
4637      return
4638      end
4639
4640
4641
4642*     **************************************
4643*     *			   	           *
4644*     *	     psi_gen_density_potentials	   *
4645*     *				           *
4646*     **************************************
4647
4648*    This routine sets the densities and potentials in psi and electron.
4649* This routine is needed for a band by band minimizer.
4650*
4651      subroutine psi_gen_density_potentials(psi_number0)
4652      implicit none
4653      integer psi_number0
4654
4655
4656#include "bafdecls.fh"
4657#include "psi.fh"
4658
4659*     ***** rhoall common block ****
4660      integer rho1_all(2)
4661      integer rho2_all(2)
4662      common / rhoall_block / rho1_all,rho2_all
4663
4664*     **** local variables ****
4665      logical zeroscf
4666      integer psi_ptr,rho_ptr,dng_ptr,rho_all_ptr,occ_ptr
4667      integer psi_number
4668
4669*     *** hacky for now ***
4670      psi_number = psi_number0
4671      zeroscf = .false.
4672      if (psi_number.gt.2) then
4673         psi_number = psi_number - 2
4674         zeroscf = .true.
4675      end if
4676      !write(*,*) "psi_number,zeroscf=",psi_number,zeroscf
4677
4678      if (psi_number.eq.1) then
4679        psi_ptr     = psi1(1)
4680        rho_ptr     = rho1(1)
4681        dng_ptr     = dng1(1)
4682        rho_all_ptr = rho1_all(1)
4683        occ_ptr     = occ1(1)
4684      else
4685        psi_ptr     = psi2(1)
4686        rho_ptr     = rho2(1)
4687        dng_ptr     = dng2(1)
4688        rho_all_ptr = rho2_all(1)
4689        occ_ptr     = occ2(1)
4690      end if
4691
4692
4693      if (zeroscf) then
4694          call electron_gen_vall0()
4695      else
4696         call electron_gen_densities(dcpl_mb(psi_ptr),
4697     >                            dbl_mb(rho_ptr),
4698     >                            dcpl_mb(dng_ptr),
4699     >                            dbl_mb(rho_all_ptr),
4700     >                            occupation_on,dbl_mb(occ_ptr))
4701         call electron_gen_scf_potentials(dbl_mb(rho_ptr),
4702     >                            dcpl_mb(dng_ptr),
4703     >                            dbl_mb(rho_all_ptr))
4704         !if (zeroscf) then
4705         !   call electron_zero_scf_potentials()
4706         !end if
4707         call electron_gen_vall()
4708      end if
4709      return
4710      end
4711
4712
4713************************ Grasmman orbitals Part ************************
4714
4715*     ***************************
4716*     *				*
4717*     *		psi_1to2	*
4718*     *				*
4719*     ***************************
4720      subroutine psi_1to2()
4721      implicit none
4722
4723#include "bafdecls.fh"
4724#include "psi.fh"
4725
4726c      call zcopy(npack1*(neq(1)+neq(2)),
4727c     >           dcpl_mb(psi1(1)),1,
4728c     >           dcpl_mb(psi2(1)),1)
4729      call Parallel_shared_vector_copy(.true.,2*npack1*(neq(1)+neq(2)),
4730     >                                 dcpl_mb(psi1(1)),
4731     >                                 dcpl_mb(psi2(1)))
4732
4733      return
4734      end
4735
4736
4737*     ***************************
4738*     *				*
4739*     *		psi_2to1	*
4740*     *				*
4741*     ***************************
4742      subroutine psi_2to1()
4743      implicit none
4744
4745#include "bafdecls.fh"
4746#include "psi.fh"
4747
4748
4749c      call zcopy(npack1*(neq(1)+neq(2)),
4750c     >           dcpl_mb(psi2(1)),1,
4751c     >           dcpl_mb(psi1(1)),1)
4752      call Parallel_shared_vector_copy(.true.,2*npack1*(neq(1)+neq(2)),
4753     >                                  dcpl_mb(psi2(1)),
4754     >                                  dcpl_mb(psi1(1)))
4755
4756c      call OrthoCheck(ispin,ne,dcpl_mb(psi1(1)))
4757      return
4758      end
4759
4760*     ***************************
4761*     *                         *
4762*     *       psiocc_1to2       *
4763*     *                         *
4764*     ***************************
4765      subroutine psiocc_1to2()
4766      implicit none
4767
4768#include "bafdecls.fh"
4769#include "psi.fh"
4770
4771      call Parallel_shared_vector_copy(.false.,2*npack1*(neq(1)+neq(2)),
4772     >                                 dcpl_mb(psi1(1)),
4773     >                                 dcpl_mb(psi2(1)))
4774      call Parallel_shared_vector_copy(.true.,(ne(1)+ne(2)),
4775     >                                 dbl_mb(occ1(1)),
4776     >                                 dbl_mb(occ2(1)))
4777      return
4778      end
4779
4780*     ***************************
4781*     *                         *
4782*     *       psiocc_2to1       *
4783*     *                         *
4784*     ***************************
4785      subroutine psiocc_2to1()
4786      implicit none
4787
4788#include "bafdecls.fh"
4789#include "psi.fh"
4790
4791      call Parallel_shared_vector_copy(.false.,2*npack1*(neq(1)+neq(2)),
4792     >                                 dcpl_mb(psi2(1)),
4793     >                                 dcpl_mb(psi1(1)))
4794      call Parallel_shared_vector_copy(.true.,(ne(1)+ne(2)),
4795     >                                 dbl_mb(occ2(1)),
4796     >                                 dbl_mb(occ1(1)))
4797      return
4798      end
4799
4800
4801*     ***************************
4802*     *                         *
4803*     *         epsi_2to1        *
4804*     *                         *
4805*     ***************************
4806      subroutine epsi_2to1()
4807      implicit none
4808
4809#include "bafdecls.fh"
4810#include "psi.fh"
4811
4812      call zcopy(npack1*(ne_excited(1)+ne_excited(2)),
4813     >           dcpl_mb(psi2_excited(1)),1,
4814     >           dcpl_mb(psi1_excited(1)),1)
4815      return
4816      end
4817
4818
4819*     ***************************
4820*     *                         *
4821*     *         epsi_1to2       *
4822*     *                         *
4823*     ***************************
4824      subroutine epsi_1to2()
4825      implicit none
4826
4827#include "bafdecls.fh"
4828#include "psi.fh"
4829
4830      call zcopy(npack1*(ne_excited(1)+ne_excited(2)),
4831     >           dcpl_mb(psi1_excited(1)),1,
4832     >           dcpl_mb(psi2_excited(1)),1)
4833      return
4834      end
4835
4836
4837
4838*     ***************************
4839*     *				*
4840*     *		psi_1get_psi	*
4841*     *				*
4842*     ***************************
4843      subroutine psi_1get_psi(rpsi)
4844      implicit none
4845      complex*16 rpsi(*)
4846
4847#include "bafdecls.fh"
4848#include "psi.fh"
4849
4850      call zcopy(npack1*(neq(1)+neq(2)),
4851     >           dcpl_mb(psi1(1)),1,
4852     >           rpsi,1)
4853
4854      return
4855      end
4856
4857
4858*     ***************************
4859*     *				*
4860*     *		psi_2get_psi	*
4861*     *				*
4862*     ***************************
4863      subroutine psi_2get_psi(rpsi)
4864      implicit none
4865      complex*16 rpsi(*)
4866
4867#include "bafdecls.fh"
4868#include "psi.fh"
4869
4870      call zcopy(npack1*(neq(1)+neq(2)),
4871     >           dcpl_mb(psi2(1)),1,
4872     >           rpsi,1)
4873
4874      return
4875      end
4876
4877*     ***************************
4878*     *				*
4879*     *		psi_check	*
4880*     *				*
4881*     ***************************
4882      subroutine psi_check()
4883      implicit none
4884
4885#include "bafdecls.fh"
4886#include "psi.fh"
4887
4888
4889      call OrthoCheck(ispin,ne,dcpl_mb(psi1(1)))
4890      return
4891      end
4892
4893
4894
4895*     ***************************
4896*     *				*
4897*     *		rho_2to1	*
4898*     *				*
4899*     ***************************
4900      subroutine rho_2to1()
4901      implicit none
4902
4903#include "bafdecls.fh"
4904#include "psi.fh"
4905
4906*     ***** rhoall common block ****
4907      integer rho1_all(2)
4908      integer rho2_all(2)
4909      common / rhoall_block / rho1_all,rho2_all
4910
4911c      call dcopy(4*nfft3d,
4912c     >           dbl_mb(rho2(1)),1,
4913c     >           dbl_mb(rho1(1)),1)
4914      call Parallel_shared_vector_copy(.true.,4*nfft3d,
4915     >                    dbl_mb(rho2(1)),
4916     >                    dbl_mb(rho1(1)))
4917
4918c      call dcopy(4*nfft3d,
4919c     >           dbl_mb(rho2_all(1)),1,
4920c     >           dbl_mb(rho1_all(1)),1)
4921      call Parallel_shared_vector_copy(.true.,4*nfft3d,
4922     >           dbl_mb(rho2_all(1)),
4923     >           dbl_mb(rho1_all(1)))
4924
4925      return
4926      end
4927
4928*     ***************************
4929*     *				*
4930*     *		rho_1to2	*
4931*     *				*
4932*     ***************************
4933      subroutine rho_1to2()
4934      implicit none
4935
4936#include "bafdecls.fh"
4937#include "psi.fh"
4938
4939*     ***** rhoall common block ****
4940      integer rho1_all(2)
4941      integer rho2_all(2)
4942      common / rhoall_block / rho1_all,rho2_all
4943
4944c      call dcopy(4*nfft3d,
4945c     >           dbl_mb(rho1(1)),1,
4946c     >           dbl_mb(rho2(1)),1)
4947c
4948c      call dcopy(4*nfft3d,
4949c     >           dbl_mb(rho1_all(1)),1,
4950c     >           dbl_mb(rho2_all(1)),1)
4951
4952      call Parallel_shared_vector_copy(.false.,4*nfft3d,
4953     >           dbl_mb(rho1(1)),
4954     >           dbl_mb(rho2(1)))
4955      call Parallel_shared_vector_copy(.true.,4*nfft3d,
4956     >           dbl_mb(rho1_all(1)),
4957     >           dbl_mb(rho2_all(1)))
4958
4959      return
4960      end
4961
4962*     ***************************
4963*     *				*
4964*     *		dng_2to1	*
4965*     *				*
4966*     ***************************
4967      subroutine dng_2to1()
4968      implicit none
4969
4970#include "bafdecls.fh"
4971#include "psi.fh"
4972
4973      call zcopy(npack0,
4974     >           dcpl_mb(dng2(1)),1,
4975     >           dcpl_mb(dng1(1)),1)
4976
4977      return
4978      end
4979
4980*     ***************************
4981*     *                         *
4982*     *         dng_1to2        *
4983*     *                         *
4984*     ***************************
4985      subroutine dng_1to2()
4986      implicit none
4987
4988#include "bafdecls.fh"
4989#include "psi.fh"
4990
4991      call zcopy(npack0,
4992     >           dcpl_mb(dng1(1)),1,
4993     >           dcpl_mb(dng2(1)),1)
4994
4995      return
4996      end
4997
4998
4999
5000*     ***********************************
5001*     *					*
5002*     *		psi_1add_oep_to_vall	*
5003*     *					*
5004*     ***********************************
5005      subroutine psi_1add_oep_to_vall()
5006      implicit none
5007
5008#include "bafdecls.fh"
5009#include "psi.fh"
5010
5011
5012      call electron_add_oep_to_vall(dbl_mb(rho1(1)))
5013
5014      return
5015      end
5016
5017
5018*     ***********************************
5019*     *					*
5020*     *		psi_1toelectron		*
5021*     *					*
5022*     ***********************************
5023      subroutine psi_1toelectron()
5024      implicit none
5025
5026#include "bafdecls.fh"
5027#include "psi.fh"
5028
5029*     ***** rhoall common block ****
5030      integer rho1_all(2)
5031      integer rho2_all(2)
5032      common / rhoall_block / rho1_all,rho2_all
5033
5034      call electron_run(dcpl_mb(psi1(1)),
5035     >                  dbl_mb(rho1(1)),
5036     >                  dcpl_mb(dng1(1)),
5037     >                  dbl_mb(rho1_all(1)),
5038     >                  occupation_on,dbl_mb(occ1(1)))
5039
5040      return
5041      end
5042
5043*     ***********************************
5044*     *					*
5045*     *		psi_1genrho		*
5046*     *					*
5047*     ***********************************
5048      subroutine psi_1genrho()
5049      implicit none
5050
5051#include "bafdecls.fh"
5052#include "psi.fh"
5053
5054*     ***** rhoall common block ****
5055      integer rho1_all(2)
5056      integer rho2_all(2)
5057      common / rhoall_block / rho1_all,rho2_all
5058
5059      call electron_genrho(dcpl_mb(psi1(1)),
5060     >                     dbl_mb(rho1(1)),
5061     >                     occupation_on,dbl_mb(occ1(1)))
5062
5063      return
5064      end
5065
5066
5067
5068*     ***********************************
5069*     *					*
5070*     *		psi_1energy		*
5071*     *					*
5072*     ***********************************
5073      real*8 function psi_1energy()
5074      implicit none
5075
5076#include "bafdecls.fh"
5077#include "psi.fh"
5078
5079*     ***** rhoall common block ****
5080      integer rho1_all(2)
5081      integer rho2_all(2)
5082      common / rhoall_block / rho1_all,rho2_all
5083
5084*     **** external functions ****
5085      real*8   electron_energy
5086      external electron_energy
5087
5088      call electron_run(dcpl_mb(psi1(1)),
5089     >                   dbl_mb(rho1(1)),
5090     >                   dcpl_mb(dng1(1)),
5091     >                   dbl_mb(rho1_all(1)),
5092     >                   occupation_on,dbl_mb(occ1(1)))
5093      psi_1energy = electron_energy(dcpl_mb(psi1(1)),
5094     >                               dbl_mb(rho1(1)),
5095     >                              dcpl_mb(dng1(1)),
5096     >                              dbl_mb(rho1_all(1)),
5097     >                              occupation_on,dbl_mb(occ1(1)))
5098
5099      return
5100      end
5101
5102*     ***********************************
5103*     *					*
5104*     *	    psi_1_noupdate_energy	*
5105*     *					*
5106*     ***********************************
5107      real*8 function psi_1_noupdate_energy()
5108      implicit none
5109
5110#include "bafdecls.fh"
5111#include "psi.fh"
5112
5113*     ***** rhoall common block ****
5114      integer rho1_all(2)
5115      integer rho2_all(2)
5116      common / rhoall_block / rho1_all,rho2_all
5117
5118*     **** external functions ****
5119      real*8   electron_energy
5120      external electron_energy
5121
5122
5123      !call electron_gen_Hpsi_k(dcpl_mb(psi1(1)))
5124      psi_1_noupdate_energy = electron_energy(dcpl_mb(psi1(1)),
5125     >                               dbl_mb(rho1(1)),
5126     >                              dcpl_mb(dng1(1)),
5127     >                              dbl_mb(rho1_all(1)),
5128     >                              occupation_on,dbl_mb(occ1(1)))
5129
5130      return
5131      end
5132
5133
5134*     ***********************************
5135*     *					*
5136*     *		psi_2energy		*
5137*     *					*
5138*     ***********************************
5139      real*8 function psi_2energy()
5140      implicit none
5141
5142#include "bafdecls.fh"
5143#include "psi.fh"
5144
5145*     ***** rhoall common block ****
5146      integer rho1_all(2)
5147      integer rho2_all(2)
5148      common / rhoall_block / rho1_all,rho2_all
5149
5150*     **** external functions ****
5151      real*8   electron_energy
5152      external electron_energy
5153
5154      call electron_run(dcpl_mb(psi2(1)),
5155     >                   dbl_mb(rho2(1)),
5156     >                  dcpl_mb(dng2(1)),
5157     >                   dbl_mb(rho2_all(1)),
5158     >                  occupation_on,dbl_mb(occ2(1)))
5159      psi_2energy = electron_energy(dcpl_mb(psi2(1)),
5160     >                               dbl_mb(rho2(1)),
5161     >                              dcpl_mb(dng2(1)),
5162     >                              dbl_mb(rho2_all(1)),
5163     >                              occupation_on,dbl_mb(occ2(1)))
5164
5165      return
5166      end
5167
5168
5169
5170*     ***********************************
5171*     *					*
5172*     *		psi_1eorbit		*
5173*     *					*
5174*     ***********************************
5175      real*8 function psi_1eorbit()
5176      implicit none
5177
5178#include "bafdecls.fh"
5179#include "psi.fh"
5180
5181*     **** external functions ****
5182      real*8   electron_eorbit
5183      external electron_eorbit
5184
5185      psi_1eorbit = electron_eorbit(dcpl_mb(psi1(1)),
5186     >                              occupation_on,dbl_mb(occ1(1)))
5187
5188      return
5189      end
5190
5191
5192*     ***********************************
5193*     *					*
5194*     *		psi_1ke 		*
5195*     *					*
5196*     ***********************************
5197      real*8 function psi_1ke()
5198      implicit none
5199
5200#include "bafdecls.fh"
5201#include "psi.fh"
5202
5203*     **** local variables ***
5204      real*8 ave
5205
5206      call ke_ave(ispin,neq,dcpl_mb(psi1(1)),ave,
5207     >            occupation_on,dbl_mb(occ1(1)))
5208
5209      psi_1ke = ave
5210      return
5211      end
5212
5213
5214
5215
5216
5217*     ***********************************
5218*     *                                 *
5219*     *         psi_1ke_atom            *
5220*     *                                 *
5221*     ***********************************
5222      real*8 function psi_1ke_atom()
5223      implicit none
5224
5225#include "bafdecls.fh"
5226#include "psi.fh"
5227
5228*     **** local variables ***
5229      real*8 ave
5230
5231*     **** external functions ****
5232      real*8   psp_kinetic_atom
5233      external psp_kinetic_atom
5234
5235      ave = psp_kinetic_atom(ispin,neq,dcpl_mb(psi1(1)))
5236
5237      psi_1ke_atom = ave
5238      return
5239      end
5240
5241*     ***********************************
5242*     *                                 *
5243*     *     psi_1valence_core_atom      *
5244*     *                                 *
5245*     ***********************************
5246      real*8 function psi_1valence_core_atom()
5247      implicit none
5248
5249#include "bafdecls.fh"
5250#include "psi.fh"
5251
5252*     **** local variables ***
5253      real*8 ave
5254
5255*     **** external functions ****
5256      real*8   psp_valence_core_atom
5257      external psp_valence_core_atom
5258
5259      ave =  psp_valence_core_atom(ispin,neq,dcpl_mb(psi1(1)))
5260
5261      psi_1valence_core_atom = ave
5262      return
5263      end
5264
5265
5266
5267*     ***********************************
5268*     *                                 *
5269*     *         psi_1vloc_atom          *
5270*     *                                 *
5271*     ***********************************
5272      real*8 function psi_1vloc_atom()
5273      implicit none
5274
5275#include "bafdecls.fh"
5276#include "psi.fh"
5277
5278*     **** local variables ***
5279      real*8 ave
5280
5281*     **** external functions ****
5282      real*8   psp_vloc_atom
5283      external psp_vloc_atom
5284
5285      ave =  psp_vloc_atom(ispin,neq,dcpl_mb(psi1(1)))
5286
5287      psi_1vloc_atom = ave
5288      return
5289      end
5290
5291
5292
5293*     ***********************************
5294*     *                                 *
5295*     *       psi_1ncmp_vloc            *
5296*     *                                 *
5297*     ***********************************
5298      real*8 function psi_1ncmp_vloc()
5299      implicit none
5300
5301#include "bafdecls.fh"
5302#include "psi.fh"
5303
5304*     **** local variables ***
5305      real*8 ave
5306
5307*     **** external functions ****
5308      real*8   psp_ncmp_vloc
5309      external psp_ncmp_vloc
5310
5311      ave =  psp_ncmp_vloc(ispin)
5312
5313      psi_1ncmp_vloc = ave
5314      return
5315      end
5316
5317
5318
5319*     ***********************************
5320*     *                                 *
5321*     *         psi_1hartree_atom      *
5322*     *                                 *
5323*     ***********************************
5324      real*8 function psi_1hartree_atom()
5325      implicit none
5326
5327#include "bafdecls.fh"
5328#include "psi.fh"
5329
5330*     **** local variables ***
5331      real*8 ave
5332
5333*     **** external functions ****
5334      real*8   psp_hartree_atom
5335      external psp_hartree_atom
5336
5337      ave =  psp_hartree_atom(ispin,neq,dcpl_mb(psi1(1)))
5338
5339      psi_1hartree_atom = ave
5340      return
5341      end
5342
5343
5344*     ***********************************
5345*     *                                 *
5346*     *      psi_1hartree_cmp_cmp       *
5347*     *                                 *
5348*     ***********************************
5349      real*8 function psi_1hartree_cmp_cmp()
5350      implicit none
5351
5352#include "psi.fh"
5353
5354*     **** external functions ****
5355      real*8   psp_hartree_cmp_cmp
5356      external psp_hartree_cmp_cmp
5357
5358      psi_1hartree_cmp_cmp = psp_hartree_cmp_cmp(ispin)
5359      return
5360      end
5361
5362*     ***********************************
5363*     *                                 *
5364*     *         psi_1hartree_cmp_pw     *
5365*     *                                 *
5366*     ***********************************
5367      real*8 function psi_1hartree_cmp_pw()
5368      implicit none
5369
5370#include "bafdecls.fh"
5371#include "psi.fh"
5372
5373*     **** external functions ****
5374      real*8   psp_hartree_cmp_pw
5375      external psp_hartree_cmp_pw
5376
5377      psi_1hartree_cmp_pw = psp_hartree_cmp_pw(ispin,dcpl_mb(dng1(1)),
5378     >                                               dbl_mb(rho1(1)))
5379      return
5380      end
5381
5382
5383c*     ***********************************
5384c*     *                                 *
5385c*     *         dng_1vlpaw_pw           *
5386c*     *                                 *
5387c*     ***********************************
5388c      real*8 function dng_1vlpaw_pw()
5389c      implicit none
5390c
5391c#include "bafdecls.fh"
5392c#include "psi.fh"
5393c
5394c*     **** external functions ****
5395c      real*8   electron_dng_vlpaw_ave
5396c      external electron_dng_vlpaw_ave
5397c
5398c      dng_1vlpaw_pw = electron_dng_vlpaw_ave(dcpl_mb(dng1(1)))
5399c
5400c      return
5401c      end
5402
5403
5404
5405*     ***********************************
5406*     *                                 *
5407*     *         psi_1xc_atom            *
5408*     *                                 *
5409*     ***********************************
5410      subroutine psi_1xc_atom(exc,pxc)
5411      implicit none
5412      real*8 exc,pxc
5413
5414#include "bafdecls.fh"
5415#include "psi.fh"
5416
5417      call psp_xc_atom(ispin,neq,dcpl_mb(psi1(1)),exc,pxc)
5418      !call D1dB_SumAll(exc)
5419      !call D1dB_SumAll(pxc)
5420      return
5421      end
5422
5423
5424*     ***********************************
5425*     *                                 *
5426*     *         psi_1qlm_atom           *
5427*     *                                 *
5428*     ***********************************
5429      subroutine psi_1qlm_atom()
5430      implicit none
5431
5432#include "bafdecls.fh"
5433#include "psi.fh"
5434
5435      call  psp_qlm_atom(ispin,neq,dcpl_mb(psi1(1)))
5436      return
5437      end
5438
5439
5440
5441
5442*     ***********************************
5443*     *					*
5444*     *		psi_1vl 		*
5445*     *					*
5446*     ***********************************
5447      real*8 function psi_1vl()
5448      implicit none
5449
5450#include "bafdecls.fh"
5451#include "psi.fh"
5452
5453
5454*     **** external functions ****
5455      real*8   electron_psi_vl_ave
5456      external electron_psi_vl_ave
5457
5458      psi_1vl = electron_psi_vl_ave(dcpl_mb(psi1(1)),dbl_mb(rho1(1)))
5459
5460      return
5461      end
5462
5463
5464
5465*     ***********************************
5466*     *                                 *
5467*     *         dng_1apc                *
5468*     *                                 *
5469*     ***********************************
5470      subroutine dng_1apc(Eapc,Papc)
5471      implicit none
5472      real*8 Eapc,Papc
5473
5474#include "bafdecls.fh"
5475#include "errquit.fh"
5476#include "psi.fh"
5477
5478*     **** local variables ****
5479      integer vtmp(2)
5480      real*8  ftmp(3)
5481
5482      if (.not.BA_push_get(mt_dcpl,npack0,'vtmp',vtmp(2),vtmp(1)))
5483     >   call errquit('dng_1apc: out of stack memory',0,MA_ERR)
5484
5485      call pspw_V_APC(ispin,ne,dcpl_mb(dng1(1)),
5486     >                         dcpl_mb(vtmp(1)),
5487     >                         Eapc,Papc,.false.,ftmp)
5488
5489      if (.not.BA_pop_stack(vtmp(2)))
5490     >  call errquit('dng_1apc popping stack',1,MA_ERR)
5491
5492      return
5493      end
5494
5495
5496*     ***********************************
5497*     *                                 *
5498*     *         dng_1vl_mm              *
5499*     *                                 *
5500*     ***********************************
5501      real*8 function dng_1vl_mm()
5502      implicit none
5503
5504#include "bafdecls.fh"
5505#include "errquit.fh"
5506#include "psi.fh"
5507
5508*     **** local variables ****
5509      integer nx,ny,nz,n2ft3d,vltmp(2),r_grid(2)
5510      real*8  elocal,esum,dv
5511
5512*     **** external functions ****
5513      integer  control_version
5514      external control_version
5515      real*8   lattice_omega
5516      external lattice_omega
5517
5518      call D3dB_n2ft3d(1,n2ft3d)
5519      call D3dB_nx(1,nx)
5520      call D3dB_ny(1,ny)
5521      call D3dB_nz(1,nz)
5522      dv = lattice_omega()/dble(nx*ny*nz)
5523
5524      if (.not.BA_push_get(mt_dbl,(n2ft3d),'vltmp',vltmp(2),vltmp(1)))
5525     >   call errquit('dng_1_vl_mm: out of stack memory',0,MA_ERR)
5526
5527*     **** average Kohn-Sham v_local energy ****
5528      call v_local_mm(dbl_mb(vltmp(1)))
5529      call Pack_cc_dot(0,dcpl_mb(dng1(1)),dbl_mb(vltmp(1)),elocal)
5530
5531*     *** add in long range part ****
5532      if (control_version().eq.4) then
5533         if (.not.BA_push_get(mt_dbl,(3*n2ft3d),'r_grid',
5534     >                       r_grid(2),r_grid(1)))
5535     >      call errquit('dng_1_vl_mm: out of stack memory',1,MA_ERR)
5536         call lattice_r_grid(dbl_mb(r_grid(1)))
5537
5538         call v_lr_local_mm(dbl_mb(r_grid(1)),dbl_mb(vltmp(1)))
5539         call D3dB_rr_dot(1,dbl_mb(rho1(1)),dbl_mb(vltmp(1)),esum)
5540         elocal = elocal + esum*dv
5541         if (.not.BA_pop_stack(r_grid(2)))
5542     >      call errquit('dng_1_vl_mm: popping stack',0,MA_ERR)
5543      end if
5544
5545      if (.not.BA_pop_stack(vltmp(2)))
5546     >  call errquit('dng_1_vl_mm: popping stack',1,MA_ERR)
5547
5548
5549      dng_1vl_mm = elocal
5550      return
5551      end
5552
5553
5554
5555
5556
5557*     ***********************************
5558*     *					*
5559*     *		psi_1vnl 		*
5560*     *					*
5561*     ***********************************
5562      real*8 function psi_1vnl()
5563      implicit none
5564
5565#include "bafdecls.fh"
5566#include "psi.fh"
5567
5568
5569*     **** external functions ****
5570      real*8   electron_psi_vnl_ave
5571      external electron_psi_vnl_ave
5572
5573      psi_1vnl = electron_psi_vnl_ave(dcpl_mb(psi1(1)),
5574     >                   occupation_on,dbl_mb(occ1(1)))
5575
5576      return
5577      end
5578
5579*     *******************************
5580*     *				    *
5581*     *		psi_1v_field 	    *
5582*     *				    *
5583*     *******************************
5584      real*8 function psi_1v_field()
5585      implicit none
5586
5587#include "bafdecls.fh"
5588#include "psi.fh"
5589
5590
5591
5592*     **** external functions ****
5593      real*8   electron_psi_v_field_ave
5594      external electron_psi_v_field_ave
5595
5596      psi_1v_field = electron_psi_v_field_ave(dcpl_mb(psi1(1)),
5597     >                                        dbl_mb(rho1(1)))
5598
5599      return
5600      end
5601
5602
5603*     ***********************************
5604*     *					*
5605*     *		rho_1Fcharge		*
5606*     *					*
5607*     ***********************************
5608      subroutine rho_1Fcharge(Fcharge)
5609      implicit none
5610      real*8 Fcharge(*)
5611
5612#include "bafdecls.fh"
5613#include "psi.fh"
5614#include "errquit.fh"
5615
5616*     **** local variables ****
5617      logical value
5618      integer n2ft3d,nx,ny,nz
5619      integer r_grid(2),rho(2)
5620      real*8  dv
5621
5622*     **** external functions ****
5623      real*8   lattice_omega
5624      external lattice_omega
5625
5626*     **** Initializationsr ****
5627      call D3dB_n2ft3d(1,n2ft3d)
5628      call D3dB_nx(1,nx)
5629      call D3dB_ny(1,ny)
5630      call D3dB_nz(1,nz)
5631      dv = lattice_omega()/dble(nx*ny*nz)
5632
5633*     **** Push memory ****
5634      value = BA_push_get(mt_dbl,(3*n2ft3d),'r_grid',
5635     >                       r_grid(2),r_grid(1))
5636      value = value.and.
5637     >        BA_push_get(mt_dbl,(3*n2ft3d),'rho',
5638     >                       rho(2),rho(1))
5639      if (.not. value) call errquit(
5640     >     'rho_1Fcharge: out of stack memory',0, MA_ERR)
5641      call Parallel_shared_vector_zero(.true.,3*n2ft3d,dbl_mb(rho(1)))
5642
5643
5644*     **** Get r_grid and rho ****
5645      call lattice_r_grid(dbl_mb(r_grid(1)))
5646      call D3dB_rr_Sum(1,dbl_mb(rho1(1)),
5647     >                   dbl_mb(rho1(1)+(ispin-1)*n2ft3d),
5648     >                   dbl_mb(rho(1)))
5649
5650*     **** Now calculate Fcharge ****
5651      call pspw_charge_rho_Fcharge(n2ft3d,dbl_mb(r_grid(1)),
5652     >                            dbl_mb(rho(1)),
5653     >                            dv,Fcharge)
5654
5655*     **** Pop memory ****
5656      value =           BA_pop_stack(rho(2))
5657      value = value.and.BA_pop_stack(r_grid(2))
5658      if (.not. value) call errquit(
5659     >     'rho_1Fcharge: error popping stack memory',0, MA_ERR)
5660
5661      return
5662      end
5663
5664
5665
5666*     ***********************************
5667*     *					*
5668*     *		rho_1exc		*
5669*     *					*
5670*     ***********************************
5671      real*8 function rho_1exc()
5672      implicit none
5673
5674#include "bafdecls.fh"
5675#include "psi.fh"
5676
5677*     ***** rhoall common block ****
5678      integer rho1_all(2)
5679      integer rho2_all(2)
5680      common / rhoall_block / rho1_all,rho2_all
5681
5682*     **** external functions ****
5683      real*8   electron_exc
5684      external electron_exc
5685
5686      rho_1exc = electron_exc(dbl_mb(rho1_all(1)))
5687      return
5688      end
5689
5690*     ***********************************
5691*     *					*
5692*     *		rho_1pxc		*
5693*     *					*
5694*     ***********************************
5695      real*8 function rho_1pxc()
5696      implicit none
5697
5698#include "bafdecls.fh"
5699#include "psi.fh"
5700
5701*     **** external functions ****
5702      real*8   electron_pxc
5703      external electron_pxc
5704
5705      rho_1pxc = electron_pxc(dbl_mb(rho1(1)))
5706      return
5707      end
5708
5709
5710*     ***********************************
5711*     *                                 *
5712*     *         psi_1meta_gga_pxc       *
5713*     *                                 *
5714*     ***********************************
5715      real*8 function psi_1meta_gga_pxc()
5716      implicit none
5717
5718#include "bafdecls.fh"
5719#include "psi.fh"
5720
5721*     **** external functions ****
5722      real*8   nwpw_meta_gga_pxc
5723      external nwpw_meta_gga_pxc
5724
5725      psi_1meta_gga_pxc = nwpw_meta_gga_pxc(ispin,neq,dcpl_mb(psi1(1)))
5726      return
5727      end
5728
5729
5730
5731*     ***********************************
5732*     *					*
5733*     *		dng_1ehartree           *
5734*     *					*
5735*     ***********************************
5736      real*8 function dng_1ehartree()
5737      implicit none
5738
5739#include "bafdecls.fh"
5740#include "psi.fh"
5741
5742*     **** external functions ****
5743      integer  control_version
5744      real*8   electron_ehartree,electron_ehartree2
5745      external control_version
5746      external electron_ehartree,electron_ehartree2
5747
5748*     **** local variables *****
5749      real*8 eh
5750
5751      eh = 0.0d0
5752      if (control_version().eq.3)
5753     >    eh = electron_ehartree(dcpl_mb(dng1(1)))
5754
5755      if (control_version().eq.4)
5756     >    eh = electron_ehartree2(dbl_mb(rho1(1)))
5757
5758      dng_1ehartree = eh
5759      return
5760      end
5761
5762*     ***********************************
5763*     *					*
5764*     * 	 psi_1vl_cosmo          *
5765*     *					*
5766*     ***********************************
5767      real*8 function psi_1vl_cosmo()
5768      implicit none
5769
5770#include "bafdecls.fh"
5771#include "psi.fh"
5772#include "errquit.fh"
5773
5774*     **** local variables ****
5775      logical ok
5776      integer r_grid(2),tmp1(2)
5777      integer nx,ny,nz,n2ft3d
5778      real*8  elocal,e1,e2,dv
5779
5780*     **** external functions ****
5781      logical  nwpw_cosmo1_on,nwpw_cosmo2_on
5782      external nwpw_cosmo1_on,nwpw_cosmo2_on
5783      real*8   lattice_omega,nwpw_cosmo_EQelcq
5784      external lattice_omega,nwpw_cosmo_EQelcq
5785
5786      if (nwpw_cosmo1_on()) then
5787         call D3dB_n2ft3d(1,n2ft3d)
5788         call D3dB_nx(1,nx)
5789         call D3dB_ny(1,ny)
5790         call D3dB_nz(1,nz)
5791         dv = lattice_omega()/dble(nx*ny*nz)
5792
5793         ok=BA_push_get(mt_dbl,(3*n2ft3d),'r_grid',r_grid(2),r_grid(1))
5794         ok=ok.and.BA_push_get(mt_dbl,(n2ft3d),'tmp1',tmp1(2),tmp1(1))
5795         if (.not.ok) call errquit("psi_1vl_cosmo: push stack",0,MA_ERR)
5796
5797         call lattice_r_grid(dbl_mb(r_grid(1)))
5798
5799         call v_local_cosmo(dbl_mb(tmp1(1)))
5800         call Pack_cc_dot(0,dcpl_mb(dng1(1)),dbl_mb(tmp1(1)),elocal)
5801
5802         call v_lr_local_cosmo(dbl_mb(r_grid(1)),dbl_mb(tmp1(1)))
5803         call D3dB_rr_dot(1,dbl_mb(rho1(1)),
5804     >                    dbl_mb(tmp1(1)),e1)
5805         call D3dB_rr_dot(1,dbl_mb(rho1(1)+(ispin-1)*n2ft3d),
5806     >                    dbl_mb(tmp1(1)),e2)
5807         elocal = elocal + (e1+e2)*dv
5808
5809         ok =        BA_pop_stack(tmp1(2))
5810         ok = ok.and.BA_pop_stack(r_grid(2))
5811         if (.not.ok) call errquit("psi_1vl_cosmo: pop stack",0,MA_ERR)
5812      else if (nwpw_cosmo2_on()) then
5813         elocal = nwpw_cosmo_EQelcq()
5814      else
5815         elocal = 0.0d0
5816      end if
5817
5818      psi_1vl_cosmo = elocal
5819      return
5820      end
5821
5822
5823
5824*     ***********************************
5825*     *					*
5826*     *		psi_2toelectron		*
5827*     *					*
5828*     ***********************************
5829      subroutine psi_2toelectron()
5830      implicit none
5831
5832#include "bafdecls.fh"
5833#include "psi.fh"
5834
5835*     ***** rhoall common block ****
5836      integer rho1_all(2)
5837      integer rho2_all(2)
5838      common / rhoall_block / rho1_all,rho2_all
5839
5840      call electron_run(dcpl_mb(psi2(1)),
5841     >                   dbl_mb(rho2(1)),
5842     >                   dcpl_mb(dng2(1)),
5843     >                   dbl_mb(rho2_all(1)),
5844     >                   occupation_on,dbl_mb(occ2(1)))
5845
5846      return
5847      end
5848
5849
5850*     ***********************************
5851*     *					*
5852*     *		psi_1check_Tangent      *
5853*     *					*
5854*     ***********************************
5855*
5856*   This routine checks the accuracy of the tangent vector.
5857*   MM = Yt*H = Yt*(I-Y*Yt)*G = Yt*G - Yt*Y*Yt*G = Yt*G - Yt*G == 0
5858
5859*     Updated - 5-18-2002
5860*
5861      subroutine psi_1check_Tangent(H)
5862      implicit none
5863      complex*16 H(*)
5864
5865#include "errquit.fh"
5866#include "bafdecls.fh"
5867#include "psi.fh"
5868
5869*     **** local variables ****
5870      logical value
5871      integer ms,n,indx,i,j
5872      integer MM(2)
5873      real*8 sum
5874
5875      do ms=1,ispin
5876         n = ne(ms)
5877         if (n.eq.0) go to 101  !*** ferromagnetic check ***
5878         value = BA_push_get(mt_dbl,n*n,'MM',MM(2),MM(1))
5879         if (.not. value)
5880     >   call errquit('out of stack memory in psi_1check_Tangent',0,
5881     &       MA_ERR)
5882
5883         indx = (ms-1)*ne(1)*npack1
5884
5885*        **** calculate MM = Yt*H ****
5886         call Grsm_ggm_dot(npack1,n,
5887     >                     dcpl_mb(psi1(1)+indx),
5888     >                     H(1+indx),
5889     >                     dbl_mb(MM(1)))
5890
5891*        **** write out MM matrix  ****
5892         sum = 0.0d0
5893         do j=1,n
5894         do i=1,n
5895            sum = sum + dbl_mb(MM(1)+(i-1)+(j-1)*n)
5896         end do
5897         end do
5898         write(*,*) "psi_1check_Tangent:",sum
5899
5900
5901
5902         value = BA_pop_stack(MM(2))
5903         if (.not. value)
5904     >    call errquit(
5905     >         'error popping stack memory in psi_1check_Tangent',0,
5906     >        MA_ERR)
5907
5908 101     continue
5909      end do
5910
5911      return
5912      end
5913
5914
5915
5916*     ***********************************
5917*     *                                 *
5918*     *         psi_2check_Tangent      *
5919*     *                                 *
5920*     ***********************************
5921*
5922*   This routine checks the accuracy of the tangent vector.
5923*   MM = Yt*H = Yt*(I-Y*Yt)*G = Yt*G - Yt*Y*Yt*G = Yt*G - Yt*G == 0
5924
5925*     Updated - 5-18-2002
5926*
5927      subroutine psi_2check_Tangent(H)
5928      implicit none
5929      complex*16 H(*)
5930
5931#include "bafdecls.fh"
5932#include "psi.fh"
5933#include "errquit.fh"
5934
5935*     **** local variables ****
5936      logical value
5937      integer ms,n,indx,i,j
5938      integer MM(2)
5939      real*8 sum
5940
5941      do ms=1,ispin
5942         n = ne(ms)
5943         if (n.eq.0) go to 101  !*** ferromagnetic check ***
5944         value = BA_push_get(mt_dbl,n*n,'MM',MM(2),MM(1))
5945         if (.not. value)
5946     >   call errquit('out of stack memory in psi_1check_Tangent',0,
5947     >        MA_ERR)
5948
5949         indx = (ms-1)*ne(1)*npack1
5950
5951*        **** calculate MM = Yt*H ****
5952         call Grsm_ggm_dot(npack1,n,
5953     >                     dcpl_mb(psi2(1)+indx),
5954     >                     H(1+indx),
5955     >                     dbl_mb(MM(1)))
5956
5957*        **** write out MM matrix  ****
5958         sum = 0.0d0
5959         do j=1,n
5960         do i=1,n
5961            sum = sum + dbl_mb(MM(1)+(i-1)+(j-1)*n)
5962         end do
5963         end do
5964         write(*,*) "psi_2check_Tangent:",sum
5965
5966
5967
5968         value = BA_pop_stack(MM(2))
5969         if (.not. value)
5970     >    call errquit(
5971     >         'error popping stack memory in psi_2check_Tangent',0,
5972     &       MA_ERR)
5973
5974 101     continue
5975      end do
5976
5977      return
5978      end
5979
5980
5981
5982*     ***********************************
5983*     *					*
5984*     *		psi_1get_Tgradient	*
5985*     *					*
5986*     ***********************************
5987
5988*     THpsi = Hpsi - Y*Y^t*Hpsi ! used by Grassman minimizers
5989*
5990      subroutine psi_1get_Tgradient(THpsi,Eout)
5991      implicit none
5992      complex*16 THpsi(*)
5993      real*8 Eout
5994
5995#include "bafdecls.fh"
5996#include "errquit.fh"
5997#include "psi.fh"
5998
5999*     ***** rhoall common block ****
6000      integer rho1_all(2)
6001      integer rho2_all(2)
6002      common / rhoall_block / rho1_all,rho2_all
6003
6004*     **** local variables ****
6005      integer tmp1(2),i,n
6006      logical value
6007
6008*     **** external functions ****
6009      logical  Dneall_m_push_get,Dneall_m_pop_stack
6010      external Dneall_m_push_get,Dneall_m_pop_stack
6011
6012      real*8   electron_energy
6013      external electron_energy
6014
6015      if (.not.Dneall_m_push_get(0,tmp1))
6016     >   call errquit('out of stack memory in psi_1get_Tradient',0,
6017     >       MA_ERR)
6018
6019      call electron_run(dcpl_mb(psi1(1)),
6020     >                   dbl_mb(rho1(1)),
6021     >                  dcpl_mb(dng1(1)),
6022     >                   dbl_mb(rho1_all(1)),
6023     >                  occupation_on,dbl_mb(occ1(1)))
6024      Eout =  electron_energy(dcpl_mb(psi1(1)),
6025     >                        dbl_mb(rho1(1)),
6026     >                        dcpl_mb(dng1(1)),
6027     >                        dbl_mb(rho1_all(1)),
6028     >                        occupation_on,dbl_mb(occ1(1)))
6029      call electron_gen_hml(dcpl_mb(psi1(1)),
6030     >                       dbl_mb(tmp1(1)))
6031
6032*     **** get the tangent THpsi = Hpsi - psi1*tmp1 , TY = HY - Y*Y'HY ****
6033      call electron_get_Tgradient(dcpl_mb(psi1(1)),
6034     >                             dbl_mb(tmp1(1)),
6035     >                            THpsi)
6036
6037      if (.not.Dneall_m_pop_stack(tmp1))
6038     > call errquit('psi_1get_Tgradient:error popping stack',1,
6039     >     MA_ERR)
6040
6041      return
6042      end
6043
6044
6045*     ***********************************
6046*     *                                 *
6047*     *         psi_1get_Tgradient0     *
6048*     *                                 *
6049*     ***********************************
6050
6051*     THpsi = Hpsi - Y*Y^t*Hpsi ! used by Grassman minimizers
6052*
6053      subroutine psi_1get_Tgradient0(THpsi,Eout)
6054      implicit none
6055      complex*16 THpsi(*)
6056      real*8 Eout
6057
6058#include "bafdecls.fh"
6059#include "errquit.fh"
6060#include "psi.fh"
6061
6062
6063*     **** local variables ****
6064      integer tmp1(2),i,n
6065      logical value
6066
6067*     **** external functions ****
6068      logical  Dneall_m_push_get,Dneall_m_pop_stack
6069      external Dneall_m_push_get,Dneall_m_pop_stack
6070
6071      real*8   electron_energy0
6072      external electron_energy0
6073
6074      if (.not.Dneall_m_push_get(0,tmp1))
6075     >   call errquit('out of stack memory in psi_1get_Tradient0',0,
6076     >       MA_ERR)
6077
6078      call electron_run0(dcpl_mb(psi1(1)))
6079      Eout =  electron_energy0(dcpl_mb(psi1(1)))
6080
6081      call electron_gen_hml(dcpl_mb(psi1(1)),
6082     >                       dbl_mb(tmp1(1)))
6083
6084*     **** get the tangent THpsi = Hpsi - psi1*tmp1 , TY = HY - Y*Y'HY ****
6085      call electron_get_Tgradient(dcpl_mb(psi1(1)),
6086     >                             dbl_mb(tmp1(1)),
6087     >                            THpsi)
6088
6089      if (.not.Dneall_m_pop_stack(tmp1))
6090     > call errquit('psi_1get_Tgradient0:error popping stack',1,
6091     >     MA_ERR)
6092
6093      return
6094      end
6095
6096
6097*     ***********************************
6098*     *                                 *
6099*     *         psi_1get_STgradient     *
6100*     *                                 *
6101*     ***********************************
6102
6103*     THpsi = Hpsi - Y*Y^t*Hpsi ! used by Grassman minimizers
6104*
6105      subroutine psi_1get_STgradient(Rpsi,THpsi,Eout)
6106      implicit none
6107      complex*16 Rpsi(*)
6108      complex*16 THpsi(*)
6109      real*8 Eout
6110
6111#include "bafdecls.fh"
6112#include "errquit.fh"
6113#include "psi.fh"
6114
6115*     ***** rhoall common block ****
6116      integer rho1_all(2)
6117      integer rho2_all(2)
6118      common / rhoall_block / rho1_all,rho2_all
6119
6120*     **** local variables ****
6121      integer tmp1(2),i,n
6122      logical value
6123
6124*     **** external functions ****
6125      logical  Dneall_m_push_get,Dneall_m_pop_stack
6126      external Dneall_m_push_get,Dneall_m_pop_stack
6127
6128      real*8   electron_energy
6129      external electron_energy
6130
6131      if (.not.Dneall_m_push_get(0,tmp1))
6132     >   call errquit('out of stack memory in psi_1get_STradient',0,
6133     >       MA_ERR)
6134
6135      call electron_run(dcpl_mb(psi1(1)),
6136     >                   dbl_mb(rho1(1)),
6137     >                  dcpl_mb(dng1(1)),
6138     >                   dbl_mb(rho1_all(1)),
6139     >                  occupation_on,dbl_mb(occ1(1)))
6140      Eout =  electron_energy(dcpl_mb(psi1(1)),
6141     >                        dbl_mb(rho1(1)),
6142     >                        dcpl_mb(dng1(1)),
6143     >                        dbl_mb(rho1_all(1)),
6144     >                        occupation_on,dbl_mb(occ1(1)))
6145      call electron_gen_hml(dcpl_mb(psi1(1)),
6146     >                       dbl_mb(tmp1(1)))
6147
6148*     **** get the tangent Rpsi = Hpsi - Spsi1*tmp1 , TY = HY - SY*Y'HY ****
6149      call psp_overlap_S(ispin,neq,dcpl_mb(psi1(1)),dcpl_mb(spsi1(1)))
6150      call electron_get_Tgradient(dcpl_mb(spsi1(1)),
6151     >                             dbl_mb(tmp1(1)),
6152     >                            Rpsi)
6153
6154*     **** THpsi = Rpsi - psi1 * <spsi1|Rpsi> ****
6155      call Grsm_gg_Copy(npack1,neq(1)+neq(2),Rpsi,THpsi)
6156      call Dneall_ffm_sym_Multiply(0,dcpl_mb(spsi1(1)),Rpsi,npack1,
6157     >                            dbl_mb(tmp1(1)))
6158      call Dneall_fmf_Multiply(0,dcpl_mb(psi1(1)),npack1,
6159     >                         dbl_mb(tmp1(1)),
6160     >                         -1.0d0,THpsi,1.0d0)
6161
6162      if (.not.Dneall_m_pop_stack(tmp1))
6163     > call errquit('psi_1get_STgradient:error popping stack',1,
6164     >     MA_ERR)
6165
6166      return
6167      end
6168
6169
6170*     ***********************************
6171*     *                                 *
6172*     *         psi_2get_STgradient     *
6173*     *                                 *
6174*     ***********************************
6175
6176*     THpsi = Hpsi - Y*Y^t*Hpsi ! used by Grassman minimizers
6177*
6178      subroutine psi_2get_STgradient(option,Rpsi,THpsi,Eout)
6179      implicit none
6180      integer    option
6181      complex*16 Rpsi(*)
6182      complex*16 THpsi(*)
6183      real*8 Eout
6184
6185#include "bafdecls.fh"
6186#include "errquit.fh"
6187#include "psi.fh"
6188
6189*     ***** rhoall common block ****
6190      integer rho1_all(2)
6191      integer rho2_all(2)
6192      common / rhoall_block / rho1_all,rho2_all
6193
6194*     **** local variables ****
6195      integer tmp1(2),i,n
6196      logical value
6197
6198*     **** external functions ****
6199      logical  Dneall_m_push_get,Dneall_m_pop_stack
6200      external Dneall_m_push_get,Dneall_m_pop_stack
6201
6202      real*8   electron_energy
6203      external electron_energy
6204
6205      if (.not.Dneall_m_push_get(0,tmp1))
6206     >   call errquit('out of stack memory in psi_2get_STradient',0,
6207     >       MA_ERR)
6208
6209      if (option.le.1) then
6210         call electron_run(dcpl_mb(psi2(1)),
6211     >                   dbl_mb(rho2(1)),
6212     >                  dcpl_mb(dng2(1)),
6213     >                   dbl_mb(rho2_all(1)),
6214     >                  occupation_on,dbl_mb(occ2(1)))
6215      end if
6216      Eout =  electron_energy(dcpl_mb(psi2(1)),
6217     >                        dbl_mb(rho2(1)),
6218     >                        dcpl_mb(dng2(1)),
6219     >                        dbl_mb(rho2_all(1)),
6220     >                        occupation_on,dbl_mb(occ2(1)))
6221      call electron_gen_hml(dcpl_mb(psi2(1)),
6222     >                       dbl_mb(tmp1(1)))
6223
6224*     **** get the tangent Rpsi = Hpsi - Spsi2*tmp1 , TY = HY - SY*Y'HY ****
6225      call psp_overlap_S(ispin,neq,dcpl_mb(psi2(1)),dcpl_mb(spsi1(1)))
6226      call electron_get_Tgradient(dcpl_mb(spsi1(1)),
6227     >                             dbl_mb(tmp1(1)),
6228     >                            Rpsi)
6229
6230*     **** THpsi = Rpsi - psi2 * <spsi2|Rpsi> ****
6231      call Grsm_gg_Copy(npack1,neq(1)+neq(2),Rpsi,THpsi)
6232      call Dneall_ffm_sym_Multiply(0,dcpl_mb(spsi1(1)),Rpsi,npack1,
6233     >                            dbl_mb(tmp1(1)))
6234      call Dneall_fmf_Multiply(0,dcpl_mb(psi2(1)),npack1,
6235     >                         dbl_mb(tmp1(1)),
6236     >                         -1.0d0,THpsi,1.0d0)
6237
6238      if (.not.Dneall_m_pop_stack(tmp1))
6239     > call errquit('psi_2get_STgradient:error popping stack',1,
6240     >     MA_ERR)
6241
6242      return
6243      end
6244
6245
6246
6247
6248
6249*     ***********************************
6250*     *					*
6251*     *		psi_1get_Gradient	*
6252*     *					*
6253*     ***********************************
6254
6255*     THpsi = Hpsi ! used by Projected Grassman minimizers
6256*
6257      subroutine psi_1get_Gradient(THpsi,Eout)
6258      implicit none
6259      complex*16 THpsi(*)
6260      real*8 Eout
6261
6262#include "bafdecls.fh"
6263#include "psi.fh"
6264
6265*     ***** rhoall common block ****
6266      integer rho1_all(2)
6267      integer rho2_all(2)
6268      common / rhoall_block / rho1_all,rho2_all
6269
6270*     **** local variables ****
6271
6272*     **** external functions ****
6273      real*8   electron_energy
6274      external electron_energy
6275
6276
6277      call electron_run(dcpl_mb(psi1(1)),
6278     >                   dbl_mb(rho1(1)),
6279     >                  dcpl_mb(dng1(1)),
6280     >                   dbl_mb(rho1_all(1)),
6281     >                   occupation_on,dbl_mb(occ1(1)))
6282
6283      Eout =  electron_energy(dcpl_mb(psi1(1)),
6284     >                        dbl_mb(rho1(1)),
6285     >                        dcpl_mb(dng1(1)),
6286     >                        dbl_mb(rho1_all(1)),
6287     >                        occupation_on,dbl_mb(occ1(1)))
6288
6289      call electron_get_Gradient(THpsi)
6290
6291      return
6292      end
6293
6294
6295*     ***********************************
6296*     *					*
6297*     *		psi_1gen_Tangent	*
6298*     *					*
6299*     ***********************************
6300
6301*     THpsi = Hpsi - Y*Y^t*Hpsi ! used by Grassman minimizers
6302*
6303      subroutine psi_1gen_Tangent(THpsi)
6304      implicit none
6305      complex*16 THpsi(*)
6306
6307#include "bafdecls.fh"
6308#include "psi.fh"
6309#include "errquit.fh"
6310
6311*     **** local variables ****
6312      integer tmp1(2)
6313
6314*     **** external functions ****
6315      logical  Dneall_m_push_get,Dneall_m_pop_stack
6316      external Dneall_m_push_get,Dneall_m_pop_stack
6317
6318      if (.not. Dneall_m_push_get(0,tmp1))
6319     >   call errquit('psi_1gen_Tangent:out of stack memory',0, MA_ERR)
6320
6321      call electron_gen_psiTangenthml(dcpl_mb(psi1(1)),
6322     >                                THpsi,
6323     >                                dbl_mb(tmp1(1)))
6324      call electron_gen_Tangent(dcpl_mb(psi1(1)),
6325     >                          dbl_mb(tmp1(1)),
6326     >                          THpsi)
6327
6328      if (.not. Dneall_m_pop_stack(tmp1))
6329     > call errquit('error popping stack memory in psi_1get_Tradient',0,
6330     &       MA_ERR)
6331
6332      return
6333      end
6334
6335
6336
6337
6338
6339*     ***********************************
6340*     *					*
6341*     *		psi_2get_Tgradient	*
6342*     *					*
6343*     ***********************************
6344      subroutine psi_2get_Tgradient(option,THpsi,Eout)
6345      implicit none
6346      integer    option
6347      complex*16 THpsi(*)
6348      real*8     Eout
6349
6350#include "errquit.fh"
6351#include "bafdecls.fh"
6352#include "psi.fh"
6353
6354*     ***** rhoall common block ****
6355      integer rho1_all(2)
6356      integer rho2_all(2)
6357      common / rhoall_block / rho1_all,rho2_all
6358
6359*     *** local variables ****
6360      integer tmp1(2)
6361
6362
6363*     **** external functions ****
6364      logical  Dneall_m_push_get,Dneall_m_pop_stack
6365      external Dneall_m_push_get,Dneall_m_pop_stack
6366
6367      real*8   electron_energy
6368      external electron_energy
6369
6370
6371!$OMP BARRIER
6372
6373      if (.not.Dneall_m_push_get(0,tmp1))
6374     >   call errquit('out of stack memory in psi_2get_Tradient',0,
6375     >       MA_ERR)
6376
6377      if (option.le.1) then
6378        call electron_run(dcpl_mb(psi2(1)),
6379     >                     dbl_mb(rho2(1)),
6380     >                    dcpl_mb(dng2(1)),
6381     >                     dbl_mb(rho2_all(1)),
6382     >                     occupation_on,dbl_mb(occ2(1)))
6383      end if
6384
6385      Eout =  electron_energy(dcpl_mb(psi2(1)),
6386     >                        dbl_mb(rho2(1)),
6387     >                        dcpl_mb(dng2(1)),
6388     >                        dbl_mb(rho2_all(1)),
6389     >                        occupation_on,dbl_mb(occ2(1)))
6390!$OMP BARRIER
6391
6392      call electron_gen_hml(dcpl_mb(psi2(1)),
6393     >                       dbl_mb(tmp1(1)))
6394      call electron_get_Tgradient(dcpl_mb(psi2(1)),
6395     >                             dbl_mb(tmp1(1)),
6396     >                             THpsi)
6397
6398      if (.not. Dneall_m_pop_stack(tmp1))
6399     >call errquit('psi_2get_Tgradient:error popping stack',1,MA_ERR)
6400
6401      return
6402      end
6403
6404
6405*     ***********************************
6406*     *                                 *
6407*     *         psi_2get_Tgradient0     *
6408*     *                                 *
6409*     ***********************************
6410      subroutine psi_2get_Tgradient0(option,THpsi,Eout)
6411      implicit none
6412      integer    option
6413      complex*16 THpsi(*)
6414      real*8     Eout
6415
6416#include "errquit.fh"
6417#include "bafdecls.fh"
6418#include "psi.fh"
6419
6420
6421*     **** local variables ****
6422      integer tmp1(2)
6423
6424*     **** external functions ****
6425      logical  Dneall_m_push_get,Dneall_m_pop_stack
6426      external Dneall_m_push_get,Dneall_m_pop_stack
6427
6428      real*8   electron_energy0
6429      external electron_energy0
6430
6431
6432!$OMP BARRIER
6433
6434      if (.not.Dneall_m_push_get(0,tmp1))
6435     >   call errquit('out of stack memory in psi_2get_Tradient0',0,
6436     >       MA_ERR)
6437
6438      if (option.le.1) then
6439        call electron_run0(dcpl_mb(psi2(1)))
6440      end if
6441
6442      Eout =  electron_energy0(dcpl_mb(psi2(1)))
6443!$OMP BARRIER
6444
6445      call electron_gen_hml(dcpl_mb(psi2(1)),
6446     >                       dbl_mb(tmp1(1)))
6447      call electron_get_Tgradient(dcpl_mb(psi2(1)),
6448     >                             dbl_mb(tmp1(1)),
6449     >                             THpsi)
6450
6451      if (.not. Dneall_m_pop_stack(tmp1))
6452     >call errquit('psi_2get_Tgradient0:error popping stack',1,MA_ERR)
6453
6454      return
6455      end
6456
6457
6458
6459*     ***********************************
6460*     *					*
6461*     *		psi_2get_Gradient	*
6462*     *					*
6463*     ***********************************
6464      subroutine psi_2get_Gradient(option,THpsi,Eout)
6465      implicit none
6466      integer    option
6467      complex*16 THpsi(*)
6468      real*8     Eout
6469
6470#include "bafdecls.fh"
6471#include "psi.fh"
6472
6473*     ***** rhoall common block ****
6474      integer rho1_all(2)
6475      integer rho2_all(2)
6476      common / rhoall_block / rho1_all,rho2_all
6477
6478*     *** local variables ****
6479
6480*     **** external functions ****
6481      real*8   electron_energy
6482      external electron_energy
6483
6484
6485      if (option.le.1) then
6486        call electron_run(dcpl_mb(psi2(1)),
6487     >                     dbl_mb(rho2(1)),
6488     >                    dcpl_mb(dng2(1)),
6489     >                     dbl_mb(rho2_all(1)),
6490     >                     occupation_on,dbl_mb(occ2(1)))
6491      end if
6492
6493      Eout =  electron_energy(dcpl_mb(psi2(1)),
6494     >                        dbl_mb(rho2(1)),
6495     >                        dcpl_mb(dng2(1)),
6496     >                        dbl_mb(rho2_all(1)),
6497     >                        occupation_on,dbl_mb(occ2(1)))
6498
6499      call electron_get_Gradient(THpsi)
6500
6501      return
6502      end
6503
6504*     ***********************************
6505*     *					*
6506*     *		psi_2gen_Tangent	*
6507*     *					*
6508*     ***********************************
6509      subroutine psi_2gen_Tangent(THpsi)
6510      implicit none
6511      complex*16 THpsi(*)
6512
6513#include "bafdecls.fh"
6514#include "psi.fh"
6515#include "errquit.fh"
6516
6517*     *** local variables ****
6518      integer tmp1(2)
6519
6520*     **** external functions ****
6521      logical  Dneall_m_push_get,Dneall_m_pop_stack
6522      external Dneall_m_push_get,Dneall_m_pop_stack
6523
6524
6525      if (.not. Dneall_m_push_get(0,tmp1))
6526     >   call errquit('psi_2gen_Tangent: out of stack memory',0,MA_ERR)
6527
6528
6529      call electron_gen_psiTangenthml(dcpl_mb(psi2(1)),
6530     >                                THpsi,
6531     >                                dbl_mb(tmp1(1)))
6532      call electron_gen_Tangent(dcpl_mb(psi2(1)),
6533     >                          dbl_mb(tmp1(1)),
6534     >                          THpsi)
6535
6536      if (.not. Dneall_m_pop_stack(tmp1))
6537     > call errquit('error popping stack memory in psi_1get_Tradient',0,
6538     &       MA_ERR)
6539
6540      return
6541      end
6542
6543
6544
6545
6546*     ***********************************
6547*     *					*
6548*     *		psi_1get_TSgradient	*
6549*     *					*
6550*     ***********************************
6551
6552*     THpsi = Hpsi - Y*Hpsi^t*Y ! used by Stiefel minimizers
6553*
6554      subroutine psi_1get_TSgradient(THpsi,Eout)
6555      implicit none
6556      complex*16 THpsi(*)
6557      real*8 Eout
6558
6559#include "errquit.fh"
6560#include "bafdecls.fh"
6561#include "psi.fh"
6562
6563*     ***** rhoall common block ****
6564      integer rho1_all(2)
6565      integer rho2_all(2)
6566      common / rhoall_block / rho1_all,rho2_all
6567
6568*     **** local variables ****
6569      integer tmp1(2)
6570
6571*     **** external functions ****
6572      logical  Dneall_m_push_get,Dneall_m_pop_stack
6573      external Dneall_m_push_get,Dneall_m_pop_stack
6574
6575      real*8   electron_energy
6576      external electron_energy
6577
6578
6579      if (.not. Dneall_m_push_get(0,tmp1))
6580     >   call errquit('psi_1get_TSradient:pushing stack',0, MA_ERR)
6581
6582
6583      call electron_run(dcpl_mb(psi1(1)),
6584     >                   dbl_mb(rho1(1)),
6585     >                  dcpl_mb(dng1(1)),
6586     >                   dbl_mb(rho1_all(1)),
6587     >                   occupation_on,dbl_mb(occ1(1)))
6588
6589      Eout =  electron_energy(dcpl_mb(psi1(1)),
6590     >                        dbl_mb(rho1(1)),
6591     >                        dcpl_mb(dng1(1)),
6592     >                        dbl_mb(rho1_all(1)),
6593     >                        occupation_on,dbl_mb(occ1(1)))
6594
6595
6596      call electron_gen_hmlt(dcpl_mb(psi1(1)),
6597     >                       dbl_mb(tmp1(1)))
6598      call electron_get_Tgradient(dcpl_mb(psi1(1)),
6599     >                             dbl_mb(tmp1(1)),
6600     >                            THpsi)
6601
6602
6603      if (.not. Dneall_m_pop_stack(tmp1))
6604     > call errquit('psi_1get_TSgradient:popping stack',1, MA_ERR)
6605
6606      return
6607      end
6608
6609
6610*     ***********************************
6611*     *					*
6612*     *		psi_2get_TSgradient	*
6613*     *					*
6614*     ***********************************
6615
6616*     THpsi = Hpsi - Y*Hpsi^t*Y ! used by Stiefel minimizers
6617*
6618      subroutine psi_2get_TSgradient(option,THpsi,Eout)
6619      implicit none
6620      integer    option
6621      complex*16 THpsi(*)
6622      real*8     Eout
6623
6624#include "errquit.fh"
6625#include "bafdecls.fh"
6626#include "psi.fh"
6627
6628*     ***** rhoall common block ****
6629      integer rho1_all(2)
6630      integer rho2_all(2)
6631      common / rhoall_block / rho1_all,rho2_all
6632
6633*     *** local variables ****
6634      integer tmp1(2)
6635
6636*     **** external functions ****
6637      logical  Dneall_m_push_get,Dneall_m_pop_stack
6638      external Dneall_m_push_get,Dneall_m_pop_stack
6639
6640      real*8   electron_energy
6641      external electron_energy
6642
6643
6644      if (.not. Dneall_m_push_get(0,tmp1))
6645     >   call errquit('psi_2get_TSgradient:pushing stack',0, MA_ERR)
6646
6647      if (option.le.1) then
6648        call electron_run(dcpl_mb(psi2(1)),
6649     >                     dbl_mb(rho2(1)),
6650     >                    dcpl_mb(dng2(1)),
6651     >                     dbl_mb(rho2_all(1)),
6652     >                    occupation_on,dbl_mb(occ2(1)))
6653      end if
6654
6655      Eout =  electron_energy(dcpl_mb(psi2(1)),
6656     >                        dbl_mb(rho2(1)),
6657     >                        dcpl_mb(dng2(1)),
6658     >                        dbl_mb(rho2_all(1)),
6659     >                        occupation_on,dbl_mb(occ2(1)))
6660
6661      call electron_gen_hmlt(dcpl_mb(psi2(1)),
6662     >                       dbl_mb(tmp1(1)))
6663      call electron_get_Tgradient(dcpl_mb(psi2(1)),
6664     >                             dbl_mb(tmp1(1)),
6665     >                             THpsi)
6666
6667      if (.not. Dneall_m_pop_stack(tmp1))
6668     > call errquit('psi_2get_TSgradient:popping stack',1, MA_ERR)
6669
6670      return
6671      end
6672
6673
6674
6675
6676*     ***********************************
6677*     *					*
6678*     *		psi_1get_TMgradient	*
6679*     *					*
6680*     ***********************************
6681      subroutine psi_1get_TMgradient(THpsi,Eout)
6682      implicit none
6683      complex*16 THpsi(*)
6684      real*8     Eout
6685
6686#include "bafdecls.fh"
6687#include "psi.fh"
6688
6689*     ***** rhoall common block ****
6690      integer rho1_all(2)
6691      integer rho2_all(2)
6692      common / rhoall_block / rho1_all,rho2_all
6693
6694*     **** external functions ****
6695      real*8   electron_energy
6696      external electron_energy
6697
6698
6699      call electron_run(dcpl_mb(psi1(1)),
6700     >                   dbl_mb(rho1(1)),
6701     >                  dcpl_mb(dng1(1)),
6702     >                   dbl_mb(rho1_all(1)),
6703     >                  occupation_on,dbl_mb(occ1(1)))
6704
6705      Eout =  electron_energy(dcpl_mb(psi1(1)),
6706     >                        dbl_mb(rho1(1)),
6707     >                        dcpl_mb(dng1(1)),
6708     >                        dbl_mb(rho1_all(1)),
6709     >                        occupation_on,dbl_mb(occ1(1)))
6710
6711      call electron_get_TMgradient(dcpl_mb(psi1(1)),
6712     >                            THpsi)
6713
6714      return
6715      end
6716
6717
6718*     ***********************************
6719*     *                                 *
6720*     *        psi_1get_TMgradient0     *
6721*     *                                 *
6722*     ***********************************
6723      subroutine psi_1get_TMgradient0(THpsi,Eout)
6724      implicit none
6725      complex*16 THpsi(*)
6726      real*8     Eout
6727
6728#include "bafdecls.fh"
6729#include "psi.fh"
6730
6731*     **** external functions ****
6732      real*8   electron_energy0
6733      external electron_energy0
6734
6735
6736      call electron_run0(dcpl_mb(psi1(1)))
6737      Eout =  electron_energy0(dcpl_mb(psi1(1)))
6738
6739      call electron_get_TMgradient(dcpl_mb(psi1(1)),
6740     >                             THpsi)
6741
6742      return
6743      end
6744
6745
6746
6747*     ***********************************
6748*     *					*
6749*     *		psi_2get_TMgradient	*
6750*     *					*
6751*     ***********************************
6752      subroutine psi_2get_TMgradient(option,THpsi,Eout)
6753      implicit none
6754      integer    option
6755      complex*16 THpsi(*)
6756      real*8     Eout
6757
6758#include "bafdecls.fh"
6759#include "psi.fh"
6760
6761*     ***** rhoall common block ****
6762      integer rho1_all(2)
6763      integer rho2_all(2)
6764      common / rhoall_block / rho1_all,rho2_all
6765
6766*     **** external functions ****
6767      real*8   electron_energy
6768      external electron_energy
6769
6770      if (option.le.1) then
6771        call electron_run(dcpl_mb(psi2(1)),
6772     >                    dbl_mb(rho2(1)),
6773     >                    dcpl_mb(dng2(1)),
6774     >                    dbl_mb(rho2_all(1)),
6775     >                    occupation_on,dbl_mb(occ2(1)))
6776      end if
6777
6778      Eout =  electron_energy(dcpl_mb(psi2(1)),
6779     >                        dbl_mb(rho2(1)),
6780     >                        dcpl_mb(dng2(1)),
6781     >                        dbl_mb(rho2_all(1)),
6782     >                        occupation_on,dbl_mb(occ2(1)))
6783
6784      call electron_get_TMgradient(dcpl_mb(psi2(1)),
6785     >                             THpsi)
6786
6787      return
6788      end
6789
6790*     ***********************************
6791*     *                                 *
6792*     *         psi_2get_TMgradient0    *
6793*     *                                 *
6794*     ***********************************
6795      subroutine psi_2get_TMgradient0(option,THpsi,Eout)
6796      implicit none
6797      integer    option
6798      complex*16 THpsi(*)
6799      real*8     Eout
6800
6801#include "bafdecls.fh"
6802#include "psi.fh"
6803
6804*     **** external functions ****
6805      real*8   electron_energy0
6806      external electron_energy0
6807
6808      if (option.le.1) then
6809        call electron_run0(dcpl_mb(psi2(1)))
6810      end if
6811
6812      Eout =  electron_energy0(dcpl_mb(psi2(1)))
6813
6814      call electron_get_TMgradient(dcpl_mb(psi2(1)),
6815     >                             THpsi)
6816
6817      return
6818      end
6819
6820
6821
6822
6823*     ***********************************
6824*     *					*
6825*     *		psi_1ke_Precondition	*
6826*     *					*
6827*     ***********************************
6828      subroutine psi_1ke_Precondition(Hpsi)
6829      implicit none
6830      complex*16 Hpsi(*)
6831
6832#include "bafdecls.fh"
6833#include "psi.fh"
6834
6835*     **** local variables ****
6836      integer neall
6837
6838      neall = neq(1)+neq(2)
6839      call ke_Precondition(npack1,neall,
6840     >                      dcpl_mb(psi1(1)),
6841     >                      Hpsi)
6842      return
6843      end
6844
6845
6846
6847*     ***********************************
6848*     *					*
6849*     *	    psi_1geodesic_transport	*
6850*     *					*
6851*     ***********************************
6852      subroutine psi_1geodesic_transport(t,H0)
6853      implicit none
6854      real*8 t
6855      complex*16 H0(*)
6856
6857#include "bafdecls.fh"
6858#include "psi.fh"
6859
6860
6861      call geodesic_transport(t,dcpl_mb(psi1(1)),H0)
6862
6863      return
6864      end
6865
6866
6867*     ***********************************
6868*     *					*
6869*     *	    psi_1geodesic_Gtransport	*
6870*     *					*
6871*     ***********************************
6872      subroutine psi_1geodesic_Gtransport(t,G0)
6873      implicit none
6874      real*8 t
6875      complex*16 G0(*)
6876
6877#include "bafdecls.fh"
6878#include "psi.fh"
6879
6880      call geodesic_Gtransport(t,dcpl_mb(psi1(1)),G0)
6881
6882      return
6883      end
6884
6885
6886
6887*     ***********************************
6888*     *                                 *
6889*     *         psi_geodesic_energy0    *
6890*     *                                 *
6891*     ***********************************
6892*
6893*    This function follows a geodesic but without updating
6894* the densities.
6895
6896      real*8 function psi_geodesic_energy0(t)
6897      implicit none
6898      real*8 t
6899
6900#include "bafdecls.fh"
6901#include "psi.fh"
6902
6903*     ***** rhoall common block ****
6904      integer rho1_all(2)
6905      integer rho2_all(2)
6906      common / rhoall_block / rho1_all,rho2_all
6907
6908      real*8 e_new
6909*     **** external functions ****
6910      real*8   electron_energy0
6911      external electron_energy0
6912
6913      call geodesic_get(t,dcpl_mb(psi1(1)),
6914     >                    dcpl_mb(psi2(1)))
6915
6916      call electron_run0(dcpl_mb(psi2(1)))
6917
6918      e_new =  electron_energy0(dcpl_mb(psi2(1)))
6919
6920      psi_geodesic_energy0 = e_new
6921      return
6922      end
6923
6924
6925
6926*     ***********************************
6927*     *					*
6928*     *		psi_geodesic_energy 	*
6929*     *					*
6930*     ***********************************
6931      real*8 function psi_geodesic_energy(t)
6932      implicit none
6933      real*8 t
6934
6935#include "bafdecls.fh"
6936#include "psi.fh"
6937
6938*     ***** rhoall common block ****
6939      integer rho1_all(2)
6940      integer rho2_all(2)
6941      common / rhoall_block / rho1_all,rho2_all
6942
6943      real*8 e_new
6944*     **** external functions ****
6945      real*8   electron_energy
6946      external electron_energy
6947
6948      call geodesic_get(t,dcpl_mb(psi1(1)),
6949     >                    dcpl_mb(psi2(1)))
6950      call electron_run(dcpl_mb(psi2(1)),
6951     >                   dbl_mb(rho2(1)),
6952     >                  dcpl_mb(dng2(1)),
6953     >                   dbl_mb(rho2_all(1)),
6954     >                   occupation_on,dbl_mb(occ2(1)))
6955      e_new =  electron_energy(dcpl_mb(psi2(1)),
6956     >                        dbl_mb(rho2(1)),
6957     >                        dcpl_mb(dng2(1)),
6958     >                        dbl_mb(rho2_all(1)),
6959     >                        occupation_on,dbl_mb(occ2(1)))
6960
6961      psi_geodesic_energy = e_new
6962      return
6963      end
6964
6965*     ***********************************
6966*     *					*
6967*     *		psi_geodesic_denergy 	*
6968*     *					*
6969*     ***********************************
6970      real*8 function psi_geodesic_denergy(t)
6971      implicit none
6972      real*8 t
6973
6974#include "bafdecls.fh"
6975#include "psi.fh"
6976
6977*     **** external functions ****
6978      real*8   electron_eorbit_noocc
6979      external electron_eorbit_noocc
6980
6981      call geodesic_transport(t,dcpl_mb(psi1(1)),
6982     >                          dcpl_mb(psi2(1)))
6983
6984      psi_geodesic_denergy
6985     > =  2.0d0*electron_eorbit_noocc(dcpl_mb(psi2(1)))
6986
6987      return
6988      end
6989
6990*     ***********************************
6991*     *					*
6992*     *		psi_geodesic_final 	*
6993*     *					*
6994*     ***********************************
6995      subroutine psi_geodesic_final(t)
6996      implicit none
6997      real*8 t
6998
6999#include "bafdecls.fh"
7000#include "psi.fh"
7001
7002      call geodesic_get(t,dcpl_mb(psi1(1)),
7003     >                    dcpl_mb(psi2(1)))
7004      return
7005      end
7006
7007
7008
7009*     ***********************************
7010*     *					*
7011*     *	    psi_1geodesic2_start	*
7012*     *					*
7013*     ***********************************
7014      subroutine psi_1geodesic2_start(H0,max_sigma,dE0)
7015      implicit none
7016      complex*16 H0(*)
7017      real*8 max_sigma
7018      real*8 dE0
7019
7020#include "bafdecls.fh"
7021#include "psi.fh"
7022
7023      call geodesic2_start(dcpl_mb(psi1(1)),H0,max_sigma,dE0)
7024
7025      return
7026      end
7027
7028*     ***********************************
7029*     *					*
7030*     *	    psi_1geodesic2_transport	*
7031*     *					*
7032*     ***********************************
7033      subroutine psi_1geodesic2_transport(t,Hnew)
7034      implicit none
7035      real*8 t
7036      complex*16 Hnew(*)
7037
7038#include "bafdecls.fh"
7039#include "psi.fh"
7040
7041      call geodesic2_transport(t,dcpl_mb(psi1(1)),Hnew)
7042
7043      return
7044      end
7045
7046
7047
7048*     ***********************************
7049*     *					*
7050*     *		psi_geodesic2_energy 	*
7051*     *					*
7052*     ***********************************
7053      real*8 function psi_geodesic2_energy(t)
7054      implicit none
7055      real*8 t
7056
7057#include "bafdecls.fh"
7058#include "psi.fh"
7059
7060*     ***** rhoall common block ****
7061      integer rho1_all(2)
7062      integer rho2_all(2)
7063      common / rhoall_block / rho1_all,rho2_all
7064
7065      real*8 e_new
7066
7067*     **** external functions ****
7068      real*8   electron_energy
7069      external electron_energy
7070
7071
7072      call geodesic2_get(t,dcpl_mb(psi1(1)),
7073     >                    dcpl_mb(psi2(1)))
7074
7075      if (occupation_on)
7076     >   call dcopy((ne(1)+ne(2)),dbl_mb(occ1(1)),1,dbl_mb(occ2(1)),1)
7077
7078*     **** check Orthogonality of psi2 **** !debug
7079*      call OrthoCheck_geo(ispin,ne,dcpl_mb(psi2(1))) !debug
7080
7081
7082      call electron_run(dcpl_mb(psi2(1)),
7083     >                   dbl_mb(rho2(1)),
7084     >                  dcpl_mb(dng2(1)),
7085     >                   dbl_mb(rho2_all(1)),
7086     >                   occupation_on,dbl_mb(occ2(1)))
7087      e_new =  electron_energy(dcpl_mb(psi2(1)),
7088     >                        dbl_mb(rho2(1)),
7089     >                        dcpl_mb(dng2(1)),
7090     >                        dbl_mb(rho2_all(1)),
7091     >                        occupation_on,dbl_mb(occ2(1)))
7092
7093      psi_geodesic2_energy = e_new
7094      return
7095      end
7096
7097*     ***********************************
7098*     *					*
7099*     *		psi_geodesic2_denergy 	*
7100*     *					*
7101*     ***********************************
7102      real*8 function psi_geodesic2_denergy(t)
7103      implicit none
7104      real*8 t
7105
7106#include "bafdecls.fh"
7107#include "psi.fh"
7108
7109*     **** external functions ****
7110      real*8   electron_eorbit
7111      external electron_eorbit
7112
7113
7114      call geodesic2_transport(t,dcpl_mb(psi1(1)),
7115     >                           dcpl_mb(psi2(1)))
7116      if (occupation_on)
7117     >   call dcopy((ne(1)+ne(2)),dbl_mb(occ1(1)),1,dbl_mb(occ2(1)),1)
7118
7119      psi_geodesic2_denergy =  2.0d0*electron_eorbit(dcpl_mb(psi2(1)),
7120     >                                  occupation_on,dbl_mb(occ2(1)))
7121
7122      return
7123      end
7124
7125*     ***********************************
7126*     *					*
7127*     *		psi_geodesic2_final 	*
7128*     *					*
7129*     ***********************************
7130      subroutine psi_geodesic2_final(t)
7131      implicit none
7132      real*8 t
7133
7134#include "bafdecls.fh"
7135#include "psi.fh"
7136
7137      integer taskid,MASTER
7138      parameter (MASTER=0)
7139c     real*8 sum1,sum2
7140
7141      call Parallel_taskid(taskid)
7142
7143      call geodesic2_get(t,dcpl_mb(psi1(1)),
7144     >                    dcpl_mb(psi2(1)))
7145      if (occupation_on)
7146     >   call dcopy((ne(1)+ne(2)),dbl_mb(occ1(1)),1,dbl_mb(occ2(1)),1)
7147      return
7148      end
7149
7150
7151
7152*     ***********************************
7153*     *					*
7154*     *		psito2_sd_update	*
7155*     *					*
7156*     ***********************************
7157      subroutine psi1to2_sd_update(dte)
7158      implicit none
7159      real*8 dte
7160
7161#include "bafdecls.fh"
7162#include "psi.fh"
7163#include "errquit.fh"
7164
7165
7166*     ***** rhoall common block ****
7167      integer rho1_all(2)
7168      integer rho2_all(2)
7169      common / rhoall_block / rho1_all,rho2_all
7170
7171*     **** local variables ****
7172      logical value
7173      integer nemaxq,ierr
7174      integer lmd(2),tmp_L(2)
7175
7176*     **** external functions ****
7177      logical  pspw_SIC,Dneall_m_push_get,Dneall_m_push_get_block
7178      logical  Dneall_m_pop_stack
7179      external pspw_SIC,Dneall_m_push_get,Dneall_m_push_get_block
7180      external Dneall_m_pop_stack
7181
7182      call electron_run(dcpl_mb(psi1(1)),
7183     >                  dbl_mb(rho1(1)),
7184     >                  dcpl_mb(dng1(1)),
7185     >                  dbl_mb(rho1_all(1)),
7186     >                  occupation_on,dbl_mb(occ1(1)))
7187
7188*     **** do a steepest descent step ****
7189      call electron_sd_update(dcpl_mb(psi1(1)),
7190     >                        dcpl_mb(psi2(1)),
7191     >                        dte)
7192
7193*     **** lagrange multiplier corrections ****
7194      nemaxq = neq(1)+neq(2)
7195
7196*     **** allocate MA local variables ****
7197      value =           Dneall_m_push_get_block(1,8,tmp_L)
7198      value = value.and.Dneall_m_push_get(0,lmd)
7199
7200c        if (occupation_on) then
7201c        call psi_lmbda2(ispin,ne,nemaxq,npack1,
7202c     >                 dcpl_mb(psi1(1)),dcpl_mb(psi2(1)),
7203c     >                 dte,dbl_mb(occ1(1)),
7204c     >                 dbl_mb(lmd(1)),
7205c     >                 dbl_mb(tmp_L(1)),ierr)
7206
7207        if (pawexist) then
7208           call psp_overlap_S(ispin,neq,
7209     >                        dcpl_mb(psi1(1)),
7210     >                        dcpl_mb(spsi1(1)))
7211           call psi_lmbda_paw(ispin,neq,nemaxq,npack1,
7212     >                        dcpl_mb(spsi1(1)),
7213     >                        dcpl_mb(psi2(1)),
7214     >                        dte,
7215     >                        dbl_mb(lmd(1)),
7216     >                        dbl_mb(tmp_L(1)),ierr)
7217        else if (occupation_on) then
7218           call psi_lmbda2(ispin,ne,nemaxq,npack1,
7219     >                     dcpl_mb(psi1(1)),
7220     >                     dcpl_mb(psi2(1)),
7221     >                     dte,dbl_mb(occ1(1)),
7222     >                     dbl_mb(lmd(1)),
7223     >                     dbl_mb(tmp_L(1)),ierr)
7224
7225        else if (pspw_SIC()) then
7226           call psi_lmbda_sic(ispin,ne,nemaxq,npack1,
7227     >                 dcpl_mb(psi1(1)),dcpl_mb(psi2(1)),dte,
7228     >                 dbl_mb(lmd(1)),
7229     >                 dbl_mb(tmp_L(1)),ierr)
7230        else
7231           call psi_lmbda(ispin,neq,nemaxq,npack1,
7232     >                 dcpl_mb(psi1(1)),dcpl_mb(psi2(1)),dte,
7233     >                 dbl_mb(lmd(1)),
7234     >                 dbl_mb(tmp_L(1)),ierr)
7235        end if
7236
7237
7238      value = value.and.Dneall_m_pop_stack(lmd)
7239      value = value.and.Dneall_m_pop_stack(tmp_L)
7240      if (.not. value)
7241     >     call errquit(
7242     >          'psi1to2_sd_update:stack failure', 0, MA_ERR)
7243      return
7244      end
7245
7246
7247*     ***********************************
7248*     *					*
7249*     *		psi_1force              *
7250*     *					*
7251*     ***********************************
7252      subroutine psi_1force(fion)
7253      implicit none
7254      real*8 fion(3,*)
7255
7256#include "bafdecls.fh"
7257#include "errquit.fh"
7258#include "psi.fh"
7259
7260*     **** local variables ****
7261      logical value
7262      integer r_grid(2),tmp(2)
7263
7264*     **** external functions ****
7265      integer  control_version
7266      external control_version
7267      logical  psp_U_psputerm,pspw_V_APC_on
7268      external psp_U_psputerm,pspw_V_APC_on
7269      logical  Dneall_m_push_get, Dneall_m_pop_stack
7270      external Dneall_m_push_get, Dneall_m_pop_stack
7271
7272c     call electron_gen_psi_r(dcpl_mb(psi1(1)))
7273c     call electron_gen_densities(dcpl_mb(psi1(1)),
7274c    >                             dbl_mb(rho1(1)),
7275c    >                            dcpl_mb(dng1(1)))
7276
7277      call f_vlocal(dcpl_mb(dng1(1)),fion)
7278
7279      if (control_version().eq.4) then
7280          value = BA_push_get(mt_dbl,(2*nfft3d),'tmp',
7281     >                        tmp(2),tmp(1))
7282          value = value.and.
7283     >            BA_push_get(mt_dbl,(6*nfft3d),'r_grid',
7284     >                        r_grid(2),r_grid(1))
7285         if (.not. value)
7286     >      call errquit('psi_1force:out of stack memory',0, MA_ERR)
7287          !call dcopy(2*nfft3d,0.0d0,0,dbl_mb(tmp(1)),1)
7288          call Parallel_shared_vector_zero(.true.,
7289     >                                     2*nfft3d,dbl_mb(tmp(1)))
7290
7291          call D3dB_rr_Sum(1,dbl_mb(rho1(1)),
7292     >                       dbl_mb(rho1(1)+(ispin-1)*2*nfft3d),
7293     >                       dbl_mb(tmp(1)))
7294          call lattice_r_grid(dbl_mb(r_grid(1)))
7295          call grad_v_lr_local(dbl_mb(r_grid(1)),
7296     >                         dbl_mb(tmp(1)),
7297     >                         fion)
7298
7299          value = BA_pop_stack(r_grid(2))
7300          value = value.and.BA_pop_stack(tmp(2))
7301         if (.not.value)
7302     >      call errquit('psi_1force:error popping stack memory',0,
7303     &                   MA_ERR)
7304      end if
7305
7306*     *** APC force here ***
7307      if (pspw_V_APC_on()) then
7308         call pspw_force_APC(ispin,ne,dcpl_mb(dng1(1)),fion)
7309      end if
7310
7311      call f_vnonlocal(ispin,
7312     >                 neq,
7313     >                 dcpl_mb(psi1(1)),
7314     >                 fion,
7315     >                 occupation_on,dbl_mb(occ1(1)))
7316
7317      if (pawexist) then
7318         value = Dneall_m_push_get(0,tmp)
7319         if (.not.value)
7320     >      call errquit('psi_1force:out of stack memory',0,MA_ERR)
7321
7322         call psi_1toelectron()
7323         call electron_gen_hml(dcpl_mb(psi1(1)),dbl_mb(tmp(1)))
7324
7325         call psp_paw_overlap_fion(ispin,
7326     >                             dbl_mb(tmp(1)),
7327     >                             dcpl_mb(psi1(1)),
7328     >                             fion)
7329
7330         value = Dneall_m_pop_stack(tmp)
7331         if (.not.value)
7332     >      call errquit('psi_1force:error popping stack',0,MA_ERR)
7333      end if
7334
7335      if (psp_U_psputerm()) then
7336         call f_psp_U_v_nonlocal(ispin,
7337     >                 neq,
7338     >                 dcpl_mb(psi1(1)),
7339     >                 fion,
7340     >                 occupation_on,dbl_mb(occ1(1)),.true.)
7341      end if
7342
7343      return
7344      end
7345
7346*     ***********************************
7347*     *                                 *
7348*     *         psi_1force_local        *
7349*     *                                 *
7350*     ***********************************
7351      subroutine psi_1force_local(fion)
7352      implicit none
7353      real*8 fion(3,*)
7354
7355#include "bafdecls.fh"
7356#include "errquit.fh"
7357#include "psi.fh"
7358
7359*     **** local variables ****
7360      logical value
7361      integer r_grid(2),tmp(2)
7362
7363*     **** external functions ****
7364      integer  control_version
7365      external control_version
7366
7367c     call electron_gen_psi_r(dcpl_mb(psi1(1)))
7368c     call electron_gen_densities(dcpl_mb(psi1(1)),
7369c    >                             dbl_mb(rho1(1)),
7370c    >                            dcpl_mb(dng1(1)))
7371
7372      call f_vlocal(dcpl_mb(dng1(1)),fion)
7373
7374      if (control_version().eq.4) then
7375          value = BA_push_get(mt_dbl,(2*nfft3d),'tmp',
7376     >                        tmp(2),tmp(1))
7377          value = value.and.
7378     >            BA_push_get(mt_dbl,(6*nfft3d),'r_grid',
7379     >                        r_grid(2),r_grid(1))
7380         if (.not. value) call errquit('out of stack memory',0, MA_ERR)
7381          call dcopy(2*nfft3d,0.0d0,0,dbl_mb(tmp(1)),1)
7382
7383          call D3dB_rr_Sum(1,dbl_mb(rho1(1)),
7384     >                       dbl_mb(rho1(1)+(ispin-1)*2*nfft3d),
7385     >                       dbl_mb(tmp(1)))
7386          call lattice_r_grid(dbl_mb(r_grid(1)))
7387          call grad_v_lr_local(dbl_mb(r_grid(1)),
7388     >                         dbl_mb(tmp(1)),
7389     >                         fion)
7390
7391          value = BA_pop_stack(r_grid(2))
7392          value = value.and.BA_pop_stack(tmp(2))
7393         if (.not. value) call errquit('error popping stack memory',0,
7394     &       MA_ERR)
7395      end if
7396
7397      return
7398      end
7399
7400*     ***********************************
7401*     *                                 *
7402*     *         psi_1force_nonlocal     *
7403*     *                                 *
7404*     ***********************************
7405      subroutine psi_1force_nonlocal(fion)
7406      implicit none
7407      real*8 fion(3,*)
7408
7409#include "bafdecls.fh"
7410#include "errquit.fh"
7411#include "psi.fh"
7412
7413      call f_vnonlocal(ispin,
7414     >                 neq,
7415     >                 dcpl_mb(psi1(1)),
7416     >                 fion,
7417     >                 occupation_on,dbl_mb(occ1(1)))
7418      return
7419      end
7420
7421
7422*     ***********************************
7423*     *                                 *
7424*     *   psi_1force_psp_U_v_nonlocal   *
7425*     *                                 *
7426*     ***********************************
7427      subroutine psi_1force_psp_U_v_nonlocal(fion)
7428      implicit none
7429      real*8 fion(3,*)
7430
7431#include "bafdecls.fh"
7432#include "errquit.fh"
7433#include "psi.fh"
7434
7435c     **** external functions ****
7436      logical  psp_U_psputerm
7437      external psp_U_psputerm
7438
7439      if (psp_U_psputerm()) then
7440         call f_psp_U_v_nonlocal(ispin,
7441     >                 neq,
7442     >                 dcpl_mb(psi1(1)),
7443     >                 fion,
7444     >                 occupation_on,dbl_mb(occ1(1)),.true.)
7445      end if
7446
7447      return
7448      end
7449
7450
7451
7452
7453
7454*     ***********************************
7455*     *					*
7456*     *		psi_1ke_stress          *
7457*     *					*
7458*     ***********************************
7459      subroutine psi_1ke_stress(stress)
7460      implicit none
7461      real*8 stress(3,3)
7462
7463#include "bafdecls.fh"
7464#include "psi.fh"
7465
7466      call ke_euv(ispin,neq,dcpl_mb(psi1(1)),stress)
7467      return
7468      end
7469
7470*     ***********************************
7471*     *					*
7472*     *		psi_1coulomb_stress     *
7473*     *					*
7474*     ***********************************
7475      subroutine psi_1coulomb_stress(stress)
7476      implicit none
7477      real*8 stress(3,3)
7478
7479#include "bafdecls.fh"
7480#include "psi.fh"
7481
7482      call coulomb_euv(dcpl_mb(dng1(1)),stress)
7483      return
7484      end
7485
7486*     ***********************************
7487*     *					*
7488*     *		rho_1exc_stress 	*
7489*     *					*
7490*     ***********************************
7491      subroutine rho_1exc_stress(stress)
7492      implicit none
7493      real*8 stress(3,3)
7494
7495#include "bafdecls.fh"
7496#include "psi.fh"
7497
7498*     ***** rhoall common block ****
7499      integer rho1_all(2)
7500      integer rho2_all(2)
7501      common / rhoall_block / rho1_all,rho2_all
7502
7503*     ***** local variables ****
7504      integer u,v,gga
7505      real*8 exc,pxc
7506      real*8 pi,scal,hm(3,3),tstress(3,3)
7507
7508*     **** external functions ****
7509      integer  control_gga
7510      real*8   rho_1exc,rho_1pxc,lattice_unitg,lattice_omega
7511      external control_gga
7512      external rho_1exc,rho_1pxc,lattice_unitg,lattice_omega
7513
7514*     *** define hm ****
7515      pi   = 4.0d0*datan(1.0d0)
7516      scal = 1.0d0/(2.0d0*pi)
7517      do v=1,3
7518      do u=1,3
7519         hm(u,v) = scal*lattice_unitg(u,v)
7520      end do
7521      end do
7522
7523*     **** LDA part ****
7524      exc = rho_1exc()
7525      pxc = rho_1pxc()
7526      do v=1,3
7527      do u=1,3
7528         stress(u,v) = (exc-pxc)*hm(u,v)
7529      end do
7530      end do
7531      !write(*,*) "hm(1,1):",hm(1,1),1.0d0/hm(1,1)
7532      !write(*,*) "exc:",exc,pxc
7533      !write(*,*) "D:",stress(1,1)
7534
7535*     **** PBE96 GGA part ****
7536*     **** finished? 11/24/04 - still need to test ***
7537      gga = control_gga()
7538      if ((gga.ge.10).and.(gga.lt.100)) then
7539       call v_bwexc_euv(gga,2*nfft3d,ispin,dbl_mb(rho1_all(1)),
7540     >                  1.0d0,1.0d0,tstress)
7541       do v=1,3
7542       do u=1,3
7543          stress(u,v) = stress(u,v) + tstress(u,v)
7544       end do
7545       end do
7546      end if
7547
7548      if (gga.eq.110) then
7549       call v_bwexc_euv(10,2*nfft3d,ispin,dbl_mb(rho1_all(1)),
7550     >                  0.75d0,1.0d0,tstress)
7551       do v=1,3
7552       do u=1,3
7553          stress(u,v) = stress(u,v) + tstress(u,v)
7554       end do
7555       end do
7556      end if
7557
7558      if (gga.eq.112) then
7559       call v_bwexc_euv(12,2*nfft3d,ispin,dbl_mb(rho1_all(1)),
7560     >                  0.75d0,1.0d0,tstress)
7561       do v=1,3
7562       do u=1,3
7563          stress(u,v) = stress(u,v) + tstress(u,v)
7564       end do
7565       end do
7566      end if
7567
7568      return
7569      end
7570
7571*     ***********************************
7572*     *					*
7573*     *		rho_1semicore_stress 	*
7574*     *					*
7575*     ***********************************
7576      subroutine rho_1semicore_stress(stress)
7577      implicit none
7578      real*8 stress(3,3)
7579
7580#include "bafdecls.fh"
7581#include "psi.fh"
7582
7583*     **** not finished ****
7584      call semicore_euv(stress)
7585
7586      return
7587      end
7588
7589
7590
7591
7592*     ***********************************
7593*     *					*
7594*     *		dng_1vlocal_stress      *
7595*     *					*
7596*     ***********************************
7597
7598      subroutine dng_1vlocal_stress(stress)
7599      implicit none
7600      real*8 stress(3,3)
7601
7602
7603#include "bafdecls.fh"
7604#include "psi.fh"
7605
7606      call v_local_euv(dcpl_mb(dng1(1)),stress)
7607
7608      return
7609      end
7610
7611*     ***********************************
7612*     *					*
7613*     *		psi_1vnonlocal_stress   *
7614*     *					*
7615*     ***********************************
7616      subroutine psi_1vnonlocal_stress(stress)
7617      implicit none
7618      real*8 stress(3,3)
7619
7620#include "bafdecls.fh"
7621#include "psi.fh"
7622
7623
7624
7625*     ***** local variables ****
7626      integer u,v
7627      real*8 evnl
7628      real*8 pi,scal,hm(3,3)
7629
7630*     **** external functions ****
7631      real*8   psi_1vnl,lattice_unitg
7632      external psi_1vnl,lattice_unitg
7633
7634*     *** define hm ****
7635      pi   = 4.0d0*datan(1.0d0)
7636      scal = 1.0d0/(2.0d0*pi)
7637      do v=1,3
7638      do u=1,3
7639         hm(u,v) = scal*lattice_unitg(u,v)
7640      end do
7641      end do
7642
7643      call v_nonlocal_euv_2(ispin,neq,dcpl_mb(psi1(1)),stress)
7644      evnl = psi_1vnl()
7645
7646      do v=1,3
7647      do u=1,3
7648         stress(u,v) = stress(u,v) - evnl*hm(u,v)
7649      end do
7650      end do
7651
7652      return
7653      end
7654
7655
7656
7657
7658*     ***********************************
7659*     *					*
7660*     *		psi_1Orb_Analysis       *
7661*     *					*
7662*     ***********************************
7663      subroutine psi_1Orb_Analysis(iunit)
7664      implicit none
7665      integer iunit
7666
7667#include "bafdecls.fh"
7668#include "psi.fh"
7669
7670c      call Orb_Analysis(iunit,ispin,ne,dcpl_mb(psi1(1)))
7671      return
7672      end
7673
7674*     ***********************************
7675*     *					*
7676*     *		psi_1Shml 	      	*
7677*     *					*
7678*     ***********************************
7679      subroutine psi_1Shml(S0,S0hml)
7680      implicit none
7681      complex*16 S0(*)
7682      complex*16 S0hml(*)
7683
7684#include "bafdecls.fh"
7685#include "psi.fh"
7686
7687      integer ms,n,shift1,shift2
7688
7689      call electron_gen_hml(dcpl_mb(psi1(1)),dbl_mb(hml(1)))
7690      do ms=1,ispin
7691            n     = ne(ms)
7692            if (n.le.0) go to 30
7693            shift1 = 1 + (ms-1)*ne(1)*npack1
7694            shift2 =     (ms-1)*ne(1)*ne(1)
7695            call DGEMM('N','N',2*npack1,n,n,
7696     >                (1.0d0),
7697     >                S0(shift1),            2*npack1,
7698     >                dbl_mb(hml(1)+shift2), n,
7699     >                (0.0d0),
7700     >                S0hml(shift1),         2*npack1)
7701   30       continue
7702      end do
7703      return
7704      end
7705
7706
7707
7708*     ***********************************
7709*     *					*
7710*     *		psi_1gen_hml      	*
7711*     *					*
7712*     ***********************************
7713      subroutine psi_1gen_hml()
7714      implicit none
7715
7716#include "bafdecls.fh"
7717#include "psi.fh"
7718
7719
7720      call electron_gen_hml(dcpl_mb(psi1(1)),dbl_mb(hml(1)))
7721
7722      return
7723      end
7724
7725
7726
7727
7728*     ***********************************
7729*     *                                 *
7730*     *         psi_1gen_hml_g          *
7731*     *                                 *
7732*     ***********************************
7733      subroutine psi_1gen_hml_g()
7734      implicit none
7735
7736#include "bafdecls.fh"
7737#include "psi.fh"
7738
7739
7740      call electron_gen_hml_g(dcpl_mb(psi1(1)),dbl_mb(hml(1)))
7741
7742      return
7743      end
7744
7745
7746*     ***********************************
7747*     *					*
7748*     *		psi_2gen_hml      	*
7749*     *					*
7750*     ***********************************
7751      subroutine psi_2gen_hml()
7752      implicit none
7753
7754#include "bafdecls.fh"
7755#include "psi.fh"
7756
7757
7758      call electron_gen_hml(dcpl_mb(psi2(1)),dbl_mb(hml(1)))
7759
7760      return
7761      end
7762
7763*     ***********************************
7764*     *					*
7765*     *		psi_hml_value     	*
7766*     *					*
7767*     ***********************************
7768      real*8  function psi_hml_value(ms,i,j)
7769      implicit none
7770      integer ms
7771      integer i,j
7772
7773#include "bafdecls.fh"
7774#include "psi.fh"
7775
7776      psi_hml_value = dbl_mb(hml(1) + (i-1)+(j-1)*ne(ms)
7777     >                              + (ms-1)*ne(1)*ne(1))
7778      return
7779      end
7780
7781*     ***********************************
7782*     *					*
7783*     *		psi_eigenvalue    	*
7784*     *					*
7785*     ***********************************
7786      real*8  function psi_eigenvalue(ms,i)
7787      implicit none
7788      integer ms
7789      integer i
7790
7791#include "bafdecls.fh"
7792#include "psi.fh"
7793
7794      real*8 sum
7795
7796      sum = dbl_mb(eig(1)+(i-1)+(ms-1)*ne(1))
7797      psi_eigenvalue = sum
7798
7799      return
7800      end
7801
7802*     ***********************************
7803*     *                                 *
7804*     *         psi_occupation          *
7805*     *                                 *
7806*     ***********************************
7807      real*8  function psi_occupation(ms,i)
7808      implicit none
7809      integer ms
7810      integer i
7811
7812#include "bafdecls.fh"
7813#include "psi.fh"
7814
7815      if (occupation_on) then
7816         psi_occupation = dbl_mb(occ1(1)+(i-1)+(ms-1)*ne(1))
7817      else
7818         psi_occupation = 1.0d0
7819      end if
7820      return
7821      end
7822
7823*     ***********************************
7824*     *                                 *
7825*     *     psi_1reverse_occupation     *
7826*     *                                 *
7827*     ***********************************
7828      subroutine psi_1reverse_occupation()
7829      implicit none
7830
7831#include "bafdecls.fh"
7832#include "psi.fh"
7833
7834      integer ms,i,indx1,indx2
7835
7836!$OMP MASTER
7837      do ms=1,ispin
7838         indx1 = occ1(1) + ne(ms) - 1 + (ms-1)*ne(1)
7839         indx2 = occ2(1)              + (ms-1)*ne(1)
7840         do i=1,ne(ms)
7841            dbl_mb(indx2)=dbl_mb(indx1)
7842            indx1 = indx1 - 1
7843            indx2 = indx2 + 1
7844         end do
7845      end do
7846      call dcopy((ne(1)+ne(2)),dbl_mb(occ2(1)),1,dbl_mb(occ1(1)),1)
7847!$OMP END MASTER
7848!$OMP BARRIER
7849      return
7850      end
7851
7852*     ***********************************
7853*     *                                 *
7854*     *     psi_1assending_occupation   *
7855*     *                                 *
7856*     ***********************************
7857*
7858*     makes sure the occupation follows an assending order
7859*
7860      subroutine psi_1assending_occupation()
7861      implicit none
7862
7863#include "bafdecls.fh"
7864#include "psi.fh"
7865
7866      if (dbl_mb(occ1(1)).lt.dbl_mb(occ1(1)+ne(1)-1)) then
7867         call psi_1reverse_occupation()
7868      end if
7869      return
7870      end
7871
7872*     ***********************************
7873*     *                                 *
7874*     *     psi_1desending_occupation   *
7875*     *                                 *
7876*     ***********************************
7877*
7878*     makes sure the occupation follows a desending order
7879*
7880      subroutine psi_1desending_occupation()
7881      implicit none
7882
7883#include "bafdecls.fh"
7884#include "psi.fh"
7885
7886      if (dbl_mb(occ1(1)+ne(1)-1).lt.dbl_mb(occ1(1))) then
7887         call psi_1reverse_occupation()
7888      end if
7889      return
7890      end
7891
7892
7893*     ***********************************
7894*     *                                 *
7895*     *      psi_0define_occupation     *
7896*     *                                 *
7897*     ***********************************
7898      subroutine psi_0define_occupation(initial_alpha,use_hml,
7899     >                                  ispin,ne,eig,hml,occ,
7900     >                                  smeartype,smearkT,
7901     >                                  smearfermi,smearcorrection)
7902      implicit none
7903      real*8  initial_alpha
7904      logical use_hml
7905      integer ispin,ne(2)
7906      real*8  eig(*),hml(*),occ(*)
7907      integer smeartype
7908      real*8 smearkT,smearfermi(2),smearcorrection
7909
7910*     **** local variables ****
7911      integer it,itmax
7912      parameter (itmax=50)
7913
7914      integer ms,nb,n,shift1,shift2,occ1_tag,ndiff
7915      real*8 e,x,kT,f,g,alpha,pi,f0
7916      real*8 ZZ,Z(2),Zlower,Zmid,Zupper,elower,emid,eupper
7917      real*8 flower,fmid,fupper,lmbda
7918
7919*     **** external functions ****
7920      integer  control_multiplicity
7921      real*8   control_TotalCharge,ion_TotalCharge_qm
7922      real*8   psi_occ_distribution,control_fractional_alpha
7923      external control_multiplicity
7924      external control_TotalCharge,ion_TotalCharge_qm
7925      external psi_occ_distribution,control_fractional_alpha
7926
7927      ZZ  = ion_TotalCharge_qm() - control_TotalCharge()
7928
7929!$OMP MASTER
7930      if (use_hml) then
7931         do ms=1,ispin
7932            shift1 = (ms-1)*ne(1)
7933            shift2 = (ms-1)*ne(1)*ne(1)
7934            do n=1,ne(ms)
7935               eig(n+shift1) = hml(n+(n-1)*ne(ms)+shift2)
7936            end do
7937         end do
7938      end if
7939      smearfermi(1)   = 0.0d0
7940      smearfermi(2)   = 0.0d0
7941      smearcorrection = 0.0d0
7942
7943
7944         if (initial_alpha.lt.0.0d0) then
7945            alpha = control_fractional_alpha()
7946         else
7947            alpha = initial_alpha
7948         end if
7949         kT    = smearkT
7950         !ZZ  = ion_TotalCharge_qm() - control_TotalCharge()
7951
7952         if (dabs(ZZ).lt.1.0d-9) go to 98
7953
7954         if (ispin.eq.2) then
7955            ndiff = control_multiplicity() - 1
7956            Z(1) = 0.5d0*(ZZ+ndiff)
7957            Z(2) = 0.5d0*(ZZ-ndiff)
7958         else
7959            Z(1) = 0.5d0*ZZ
7960            Z(2) = 0.0d0
7961         end if
7962
7963         pi    = 4.0d0*datan(1.0d0)
7964         !if (initial) alpha = 1.0d0
7965         if (smeartype.le.0) alpha = 0.0d0
7966
7967*        **** outer loop over spins ****
7968         smearcorrection = 0.0d0
7969         do ms=1,ispin
7970
7971*           **** find eupper and elower ****
7972            elower =  9.9d12
7973            eupper = -9.9d12
7974            shift1 = (ms-1)*ne(1) + 1
7975            do n=1,ne(ms)
7976              e       = eig(shift1)
7977              if (e.lt.elower) elower = e
7978              if (e.gt.eupper) eupper = e
7979              shift1  = shift1 + 1
7980            end do
7981
7982*           **** find fermi level ****
7983            Zlower = 0.0d0
7984            Zupper = 0.0d0
7985            shift1 = (ms-1)*ne(1) + 1
7986            do n=1,ne(ms)
7987              e = eig(shift1)
7988              Zlower = Zlower
7989     >         + psi_occ_distribution(smeartype,(e-elower)/kT)
7990              Zupper = Zupper
7991     >         + psi_occ_distribution(smeartype,(e-eupper)/kT)
7992              shift1  = shift1 + 1
7993            end do
7994
7995
7996            flower = Zlower - Z(ms)
7997            fupper = Zupper - Z(ms)
7998
7999            if (flower*fupper.ge.0.0d0)
8000     >       call errquit(
8001     >            'psi_0define_occupation:Fermi energy not found',ms,0)
8002
8003            it = 0
8004  20        it = it + 1
8005            emid = 0.5d0*(elower + eupper)
8006            Zmid = 0.0d0
8007            shift1 = (ms-1)*ne(1) + 1
8008            do n=1,ne(ms)
8009              e = eig(shift1)
8010              Zmid = Zmid + psi_occ_distribution(smeartype,(e-emid)/kT)
8011              shift1  = shift1 + 1
8012            end do
8013            fmid = Zmid - Z(ms)
8014            if (fmid.lt.0.0d0) then
8015               flower = fmid
8016               elower = emid
8017            else
8018               fupper = fmid
8019               eupper = emid
8020            end if
8021            if ( (dabs(fmid)     .gt.1.0d-11) .and.
8022     >           ((eupper-elower).gt.1.0d-11) .and.
8023     >           (it.lt.itmax))   goto 20
8024
8025
8026            smearfermi(ms) = emid
8027
8028*           **** determine filling and correction ****
8029            shift1 = (ms-1)*ne(1) + 1
8030            shift2 = (ms-1)*ne(1) + 1
8031            do n=1,ne(ms)
8032              e = eig(shift1)
8033              x = (e - smearfermi(ms))/kT
8034              f = psi_occ_distribution(smeartype,x)
8035              f0 = occ(shift2)
8036
8037              occ(shift2) = (1.0d0-alpha)*f0 + alpha*f
8038
8039              if (smeartype.eq.1) then
8040                 if (  (occ(shift2)       .gt.1.0d-6) .and.
8041     >               ( (1.0d0-occ(shift2)).gt.1.0d-6)) then
8042                smearcorrection = smearcorrection
8043     >           + kT*( occ(shift2)*log(occ(shift2))
8044     >           + (1.0d0-occ(shift2))*log(1.0d0-occ(shift2)) )
8045                 end if
8046              else if (smeartype.eq.2) then
8047                smearcorrection
8048     >              = smearcorrection
8049     >              - kT*dexp(-x*x)/(4.0d0*dsqrt(pi))
8050              else if (smeartype.eq.4) then
8051                smearcorrection
8052     >              = smearcorrection
8053     >              - kT*dexp(-(x+dsqrt(0.5d0))*(x+dsqrt(0.5d0)))
8054     >              * (1.0d0 + dsqrt(2.0d0) * x)
8055     >              / (2.0d0 * dsqrt(pi))
8056              end if
8057
8058              shift1  = shift1 + 1
8059              shift2  = shift2 + 1
8060            end do
8061
8062         end do !** ms***
8063         if (ms.eq.1) smearcorrection=smearcorrection+smearcorrection
8064
8065         go to 99
8066
8067
8068  98     continue
8069         smearcorrection = 0.0d0
8070         call dcopy(ne(1)+ne(2),0.0d0,0,occ,1)
8071
8072  99     continue
8073
8074!$OMP END MASTER
8075!$OMP BARRIER
8076
8077      return
8078      end
8079
8080
8081
8082
8083
8084
8085*     ***********************************
8086*     *                                 *
8087*     *      psi_1define_occupation     *
8088*     *                                 *
8089*     ***********************************
8090      subroutine psi_1define_occupation(initial_alpha,use_hml)
8091      implicit none
8092      real*8  initial_alpha
8093      logical use_hml
8094
8095#include "bafdecls.fh"
8096#include "psi.fh"
8097
8098*     **** local variables ****
8099      integer it,itmax
8100      parameter (itmax=50)
8101
8102      integer ms,nb,n,shift1,shift2,occ1_tag,ndiff
8103      real*8 e,x,kT,f,g,alpha,pi,f0
8104      real*8 ZZ,Z(2),Zlower,Zmid,Zupper,elower,emid,eupper
8105      real*8 flower,fmid,fupper,lmbda
8106
8107*     **** external functions ****
8108      integer  control_multiplicity
8109      real*8   control_TotalCharge,ion_TotalCharge_qm
8110      real*8   psi_occ_distribution,control_fractional_alpha
8111      external control_multiplicity
8112      external control_TotalCharge,ion_TotalCharge_qm
8113      external psi_occ_distribution,control_fractional_alpha
8114
8115      ZZ  = ion_TotalCharge_qm() - control_TotalCharge()
8116
8117!$OMP MASTER
8118      if (use_hml) then
8119         do ms=1,ispin
8120            shift1 = eig(1) + (ms-1)*ne(1)
8121            shift2 = hml(1) + (ms-1)*ne(1)*ne(1)
8122            do n=1,ne(ms)
8123               dbl_mb(shift1+n-1) = dbl_mb(shift2+(n-1)+(n-1)*ne(ms))
8124            end do
8125         end do
8126      end if
8127      smearfermi(1)   = 0.0d0
8128      smearfermi(2)   = 0.0d0
8129      smearcorrection = 0.0d0
8130
8131      if (occupation_on) then
8132         if (initial_alpha.lt.0.0d0) then
8133            alpha = control_fractional_alpha()
8134         else
8135            alpha = initial_alpha
8136         end if
8137         kT    = smearkT
8138         !ZZ  = ion_TotalCharge_qm() - control_TotalCharge()
8139
8140         if (dabs(ZZ).lt.1.0d-9) go to 98
8141
8142         if (ispin.eq.2) then
8143            ndiff = control_multiplicity() - 1
8144            Z(1) = 0.5d0*(ZZ+ndiff)
8145            Z(2) = 0.5d0*(ZZ-ndiff)
8146         else
8147            Z(1) = 0.5d0*ZZ
8148            Z(2) = 0.0d0
8149         end if
8150
8151         pi    = 4.0d0*datan(1.0d0)
8152         !if (initial) alpha = 1.0d0
8153         if (smeartype.le.0) alpha = 0.0d0
8154
8155*        **** outer loop over spins ****
8156         smearcorrection = 0.0d0
8157         do ms=1,ispin
8158
8159*           **** find eupper and elower ****
8160            elower =  9.9d12
8161            eupper = -9.9d12
8162            shift1 = eig(1) + (ms-1)*ne(1)
8163            do n=1,ne(ms)
8164              e       = dbl_mb(shift1)
8165              if (e.lt.elower) elower = e
8166              if (e.gt.eupper) eupper = e
8167              shift1  = shift1 + 1
8168            end do
8169
8170
8171*           **** find fermi level ****
8172            Zlower = 0.0d0
8173            Zupper = 0.0d0
8174            shift1 = eig(1) + (ms-1)*ne(1)
8175            do n=1,ne(ms)
8176              e = dbl_mb(shift1)
8177              Zlower = Zlower
8178     >         + psi_occ_distribution(smeartype,(e-elower)/kT)
8179              Zupper = Zupper
8180     >         + psi_occ_distribution(smeartype,(e-eupper)/kT)
8181              shift1  = shift1 + 1
8182            end do
8183
8184            flower = Zlower - Z(ms)
8185            fupper = Zupper - Z(ms)
8186
8187            if (flower*fupper.ge.0.0d0)
8188     >       call errquit(
8189     >            'psi_1define_occupation:Fermi energy not found',ms,0)
8190
8191            it = 0
8192  20        it = it + 1
8193            emid = 0.5d0*(elower + eupper)
8194            Zmid = 0.0d0
8195            shift1 = eig(1) + (ms-1)*ne(1)
8196            do n=1,ne(ms)
8197              e = dbl_mb(shift1)
8198              Zmid = Zmid + psi_occ_distribution(smeartype,(e-emid)/kT)
8199              shift1  = shift1 + 1
8200            end do
8201            fmid = Zmid - Z(ms)
8202            if (fmid.lt.0.0d0) then
8203               flower = fmid
8204               elower = emid
8205            else
8206               fupper = fmid
8207               eupper = emid
8208            end if
8209            if ( (dabs(fmid)     .gt.1.0d-11) .and.
8210     >           ((eupper-elower).gt.1.0d-11) .and.
8211     >           (it.lt.itmax))   goto 20
8212
8213
8214            smearfermi(ms) = emid
8215
8216*           **** determine filling and correction ****
8217            shift1 = eig(1)  + (ms-1)*ne(1)
8218            shift2 = occ1(1) + (ms-1)*ne(1)
8219            do n=1,ne(ms)
8220              e = dbl_mb(shift1)
8221              x = (e - smearfermi(ms))/kT
8222              f = psi_occ_distribution(smeartype,x)
8223              f0 = dbl_mb(shift2)
8224
8225              dbl_mb(shift2) = (1.0d0-alpha)*f0 + alpha*f
8226
8227              if (smeartype.eq.1) then
8228                 if (  (dbl_mb(shift2)       .gt.1.0d-6) .and.
8229     >               ( (1.0d0-dbl_mb(shift2)).gt.1.0d-6)) then
8230                smearcorrection = smearcorrection
8231     >           + kT*( dbl_mb(shift2)*log(dbl_mb(shift2))
8232     >           + (1.0d0-dbl_mb(shift2))*log(1.0d0-dbl_mb(shift2)) )
8233                 end if
8234              else if (smeartype.eq.2) then
8235                smearcorrection
8236     >              = smearcorrection
8237     >              - kT*dexp(-x*x)/(4.0d0*dsqrt(pi))
8238              else if (smeartype.eq.4) then
8239                smearcorrection
8240     >              = smearcorrection
8241     >              - kT*dexp(-(x+dsqrt(0.5d0))*(x+dsqrt(0.5d0)))
8242     >              * (1.0d0 + dsqrt(2.0d0) * x)
8243     >              / (2.0d0 * dsqrt(pi))
8244              end if
8245
8246              shift1  = shift1 + 1
8247              shift2  = shift2 + 1
8248            end do
8249
8250         end do !** ms***
8251         if (ms.eq.1) smearcorrection=smearcorrection+smearcorrection
8252
8253         go to 99
8254
8255  98     continue
8256         smearcorrection = 0.0d0
8257         do ms=1,ispin
8258            shift2 = occ1(1) + (ms-1)*ne(1)
8259            call dcopy(ne(ms),0.0d0,0,dbl_mb(shift2),1)
8260         end do
8261
8262  99     continue
8263
8264      end if
8265!$OMP END MASTER
8266!$OMP BARRIER
8267
8268      return
8269      end
8270
8271
8272c  set nwpw:fractional_smeartype 1 #0-none, 1-Fermi-Dirac, 2-Gaussian, 3-Hermite
8273c                                   4-Marzari-Vanderbilt
8274      real*8 function psi_occ_distribution(smeartype,e)
8275      implicit none
8276      integer smeartype
8277      real*8 e
8278      real*8 f
8279
8280*     **** external functions ****
8281      real*8   util_erfc
8282      external util_erfc
8283
8284      if (smeartype.eq.1) then
8285         if (e.gt.30.0d0) then
8286           f = 0.0d0
8287         else if (e.lt.(-30.0d0)) then
8288           f = 1.0d0
8289         else
8290           f = 1.0d0/(1.0d0+dexp(e))
8291         end if
8292      else if (smeartype.eq.2) then
8293         f = 0.5d0*util_erfc(e)
8294      else if (smeartype.eq.4) then
8295         f = dexp(-(e + dsqrt(0.5d0)) * (e + dsqrt(0.5d0)))
8296     >     * dsqrt(0.125d0 / datan(1.0d0))
8297     >     + 0.5d0 * util_erfc(e + dsqrt(0.5d0))
8298      else
8299         if (e.gt.0.0d0) then
8300           f = 0.0d0
8301         else
8302           f = 1.0d0
8303         end if
8304      end if
8305      psi_occ_distribution = f
8306      return
8307      end
8308
8309      real*8 function psi_smearfermi(ms)
8310      implicit none
8311      integer ms
8312#include "psi.fh"
8313      psi_smearfermi = smearfermi(ms)
8314      return
8315      end
8316      real*8 function psi_smearcorrection()
8317      implicit none
8318#include "psi.fh"
8319      psi_smearcorrection = smearcorrection
8320      return
8321      end
8322
8323
8324
8325
8326
8327*     ***********************************
8328*     *                                 *
8329*     *          psi_virtual            *
8330*     *                                 *
8331*     ***********************************
8332      real*8  function psi_virtual(ms,i)
8333      implicit none
8334      integer ms
8335      integer i
8336
8337#include "bafdecls.fh"
8338#include "psi.fh"
8339
8340      psi_virtual=dbl_mb(eig_excited(1)+(i-1)+(ms-1)*ne_excited(1))
8341
8342      return
8343      end
8344
8345*     ***********************************
8346*     *					*
8347*     *		psi_hml		   	*
8348*     *					*
8349*     ***********************************
8350      real*8  function psi_hml(ms,i,j)
8351      implicit none
8352      integer ms
8353      integer i,j
8354
8355#include "bafdecls.fh"
8356#include "psi.fh"
8357
8358      psi_hml = dbl_mb(hml(1)-1 + i
8359     >                          + (j-1)*ne(ms)
8360     >                          + (ms-1)*ne(1)*ne(1))
8361
8362      return
8363      end
8364
8365
8366*     ***********************************
8367*     *                                 *
8368*     *         psi_iptr_hml            *
8369*     *                                 *
8370*     ***********************************
8371      integer function psi_iptr_hml(ms,i,j)
8372      implicit none
8373      integer ms
8374      integer i,j
8375
8376#include "bafdecls.fh"
8377#include "psi.fh"
8378
8379      psi_iptr_hml = (hml(1)-1 + i
8380     >                          + (j-1)*ne(ms)
8381     >                          + (ms-1)*ne(1)*ne(1))
8382
8383      return
8384      end
8385
8386
8387*     ***********************************
8388*     *					*
8389*     *		psi_spin_density  	*
8390*     *					*
8391*     ***********************************
8392      subroutine psi_spin_density(en)
8393      implicit none
8394      real*8 en(2)
8395
8396#include "bafdecls.fh"
8397#include "psi.fh"
8398
8399*     **** local variables ****
8400      integer ms,nx,ny,nz,n2ft3d
8401      real*8  scale,sumall
8402
8403*     **** external functions ****
8404      real*8   lattice_omega
8405      external lattice_omega
8406
8407      call D3dB_n2ft3d(1,n2ft3d)
8408      !n2ft3d = 2*n2ft3d
8409      call D3dB_nx(1,nx)
8410      call D3dB_ny(1,ny)
8411      call D3dB_nz(1,nz)
8412      scale = lattice_omega()/dble(nx*ny*nz)
8413
8414*     **** check total number of electrons ****
8415      en(2) = 0.0d0
8416      do ms =1,ispin
8417         call D3dB_r_dsum(1,dbl_mb(rho1(1)+(ms-1)*n2ft3d),sumall)
8418         en(ms) = sumall*scale
8419      end do
8420
8421      return
8422      end
8423
8424*     ***********************************
8425*     *					*
8426*     *		psi_spin2     	        *
8427*     *					*
8428*     ***********************************
8429      subroutine psi_spin2(Sab)
8430      implicit none
8431      real*8 Sab
8432
8433#include "bafdecls.fh"
8434#include "psi.fh"
8435
8436      call Calculate_psi_spin2(ispin,ne,npack1,dcpl_mb(psi1(1)),
8437     >                         occupation_on,dbl_mb(occ1(1)),Sab)
8438      return
8439      end
8440
8441*     ***********************************
8442*     *					*
8443*     *		psi_1rotate2       	*
8444*     *					*
8445*     ***********************************
8446      subroutine psi_1rotate2()
8447      implicit none
8448
8449#include "bafdecls.fh"
8450#include "psi.fh"
8451
8452c*     ***** local variables *****
8453c      integer ms,index,i,j,shift1,shift2
8454
8455      call Dneall_fmf_Multiply(0,dcpl_mb(psi1(1)),npack1,
8456     >                          dbl_mb(hml(1)),1.0d0,
8457     >                          dcpl_mb(psi2(1)),0.0d0)
8458
8459
8460
8461c      !call dcopy(2*npack1*(ne(1)+ne(2)),0.0d0,0,dcpl_mb(psi2(1)),1)
8462c      do ms=1,ispin
8463c         if (ne(ms).le.0) go to 30
8464c         shift1 = (ms-1)*ne(1)
8465c         shift2 = (ms-1)*ne(1)*ne(1)
8466c
8467c         call DGEMM('N','N',2*npack1,ne(ms),ne(ms),
8468c     >              (1.0d0),
8469c     >              dcpl_mb(psi1(1)+shift1*npack1),2*npack1,
8470c     >              dbl_mb(hml(1)+shift2),ne(ms),
8471c     >              (0.0d0),
8472c     >              dcpl_mb(psi2(1)+shift1*npack1),2*npack1)
8473cc        do j=1,ne(ms)
8474cc          do i=1,ne(ms)
8475cc             index = (i-1) + (j-1)*ne(ms) + shift2
8476cc
8477cc              call D3dB_cc_daxpy(1,dbl_mb(hml(1)+index),
8478cc     >                           dcpl_mb(psi1(1)+(i-1+shift1)*nfft3d),
8479cc     >                           dcpl_mb(psi2(1)+(j-1+shift1)*nfft3d))
8480cc             call Pack_cc_daxpy(1,dbl_mb(hml(1)+index),
8481cc    >                           dcpl_mb(psi1(1)+(i-1+shift1)*npack1),
8482cc    >                           dcpl_mb(psi2(1)+(j-1+shift1)*npack1))
8483cc          end do
8484cc        end do
8485c
8486c   30   continue
8487c      end do
8488
8489      return
8490      end
8491
8492*     ***********************************
8493*     *					*
8494*     *		psi_2rotate1       	*
8495*     *					*
8496*     ***********************************
8497      subroutine psi_2rotate1()
8498      implicit none
8499
8500#include "bafdecls.fh"
8501#include "psi.fh"
8502
8503c*     ***** local variables *****
8504c      integer ms,index,i,j,shift1,shift2
8505
8506      call Dneall_fmf_Multiply(0,dcpl_mb(psi2(1)),npack1,
8507     >                          dbl_mb(hml(1)),1.0d0,
8508     >                          dcpl_mb(psi1(1)),0.0d0)
8509
8510c      do ms=1,ispin
8511c         if (ne(ms).le.0) go to 30
8512c         shift1 = (ms-1)*ne(1)
8513c         shift2 = (ms-1)*ne(1)*ne(1)
8514c
8515c         call DGEMM('N','N',2*npack1,ne(ms),ne(ms),
8516c     >              (1.0d0),
8517c     >              dcpl_mb(psi2(1)+shift1*npack1),2*npack1,
8518c     >              dbl_mb(hml(1)+shift2),ne(ms),
8519c     >              (0.0d0),
8520c     >              dcpl_mb(psi1(1)+shift1*npack1),2*npack1)
8521c
8522c   30    continue
8523c      end do
8524
8525      return
8526      end
8527
8528
8529*     ***********************************
8530*     *					*
8531*     *		psi_diagonalize_hml	*
8532*     *					*
8533*     ***********************************
8534      subroutine psi_diagonalize_hml()
8535      implicit none
8536
8537#include "bafdecls.fh"
8538#include "psi.fh"
8539
8540
8541c*     ***** local variables ****
8542c      logical value
8543c      integer ms,shift1,shift2,ierr,i,j,indx
8544c      integer tmp1(2)
8545
8546
8547      call Dneall_m_diagonalize(0,dbl_mb(hml(1)),
8548     >                             dbl_mb(eig(1)),.false.)
8549
8550c      value = BA_push_get(mt_dbl,(2*ne(1)*ne(1)),'tmp1',tmp1(2),tmp1(1))
8551c      if (.not. value)
8552c     >   call errquit('out of stack memory in psi_diagonalize_hml',0,
8553c     &       MA_ERR)
8554
8555
8556c*     ***** diagonalize the hamiltonian matrix *****
8557c      call dcopy((ne(1)+ne(2)),0.0d0,0,dbl_mb(eig(1)),1)
8558c      do ms=1,ispin
8559c         shift1 = (ms-1)*ne(1)
8560c         shift2 = (ms-1)*ne(1)*ne(1)
8561c         if (ne(ms).le.0) go to 30
8562
8563cc        call eigen(ne(ms),ne(ms),
8564cc    >              dbl_mb(hml(1)+shift2),
8565cc    >              dbl_mb(eig(1)+shift1),
8566cc    >              dbl_mb(tmp1(1)),ierr)
8567c
8568c         call DSYEV('V','U',ne(ms),
8569c     >              dbl_mb(hml(1)+shift2),ne(ms),
8570c     >              dbl_mb(eig(1)+shift1),
8571c     >              dbl_mb(tmp1(1)),2*ne(1)*ne(1),
8572c     >              ierr)
8573c
8574c         call eigsrt(dbl_mb(eig(1)+shift1),
8575c     >              dbl_mb(hml(1)+shift2),
8576c     >              ne(ms),ne(ms))
8577c
8578c  30    continue
8579c      end do
8580c
8581c
8582c      value = BA_pop_stack(tmp1(2))
8583c      if (.not. value)
8584c     > call errquit('error popping stack in psi_diagonalize_hml',0,
8585c     &       MA_ERR)
8586
8587      return
8588      end
8589
8590*     ***********************************
8591*     *					*
8592*     *	  psi_diagonalize_hml_assending *
8593*     *					*
8594*     ***********************************
8595      subroutine psi_diagonalize_hml_assending()
8596      implicit none
8597
8598#include "bafdecls.fh"
8599#include "psi.fh"
8600
8601
8602c*     ***** local variables ****
8603c      logical value
8604c      integer ms,shift1,shift2,ierr
8605c      integer tmp1(2)
8606
8607      call Dneall_m_diagonalize(0,dbl_mb(hml(1)),dbl_mb(eig(1)),.true.)
8608
8609c      value = BA_push_get(mt_dbl,(2*ne(1)*ne(1)),'tmp1',tmp1(2),tmp1(1))
8610c      if (.not. value)
8611c     >   call errquit(
8612c     >    'out of stack memory in psi_diagonalize_hml_assending',0,
8613c     >     MA_ERR)
8614c
8615c
8616c*     ***** diagonalize the hamiltonian matrix *****
8617c      call dcopy((ne(1)+ne(2)),0.0d0,0,dbl_mb(eig(1)),1)
8618c      do ms=1,ispin
8619c         shift1 = (ms-1)*ne(1)
8620c         shift2 = (ms-1)*ne(1)*ne(1)
8621c         if (ne(ms).le.0) go to 30
8622c
8623c         call DSYEV('V','U',ne(ms),
8624c     >              dbl_mb(hml(1)+shift2),ne(ms),
8625c     >              dbl_mb(eig(1)+shift1),
8626c     >              dbl_mb(tmp1(1)),2*ne(1)*ne(1),
8627c     >              ierr)
8628c
8629c   30    continue
8630c      end do
8631
8632
8633c      value = BA_pop_stack(tmp1(2))
8634c      if (.not. value)
8635c     > call errquit(
8636c     >   'error popping stack in psi_diagonalize_hml_assending',0,
8637c     >     MA_ERR)
8638
8639      return
8640      end
8641
8642
8643
8644*     ***************************
8645*     *				*
8646*     *		psi_error	*
8647*     *				*
8648*     ***************************
8649      real*8 function psi_error()
8650      implicit none
8651#include "errquit.fh"
8652
8653#include "bafdecls.fh"
8654#include "psi.fh"
8655
8656*     ***** local variables ****
8657      logical value
8658      integer k,n
8659      real*8  error,sum,size
8660      integer tmp1(2)
8661
8662      value = BA_push_get(mt_dcpl,(npack1),'tmp1',tmp1(2),tmp1(1))
8663      if (.not. value)
8664     >   call errquit('out of stack memory in psi_error',0, MA_ERR)
8665
8666
8667      error = 0.0d0
8668      size =  dble(ne(1)+ne(2))
8669      do n=1, (neq(1)+neq(2))
8670         do k=1,npack1
8671            dcpl_mb(tmp1(1)+k-1) = dcpl_mb(psi2(1)+k-1+(n-1)*npack1)
8672     >                           - dcpl_mb(psi1(1)+k-1+(n-1)*npack1)
8673         end do
8674c         call D3dB_cc_dot(1,dcpl_mb(tmp1(1)),dcpl_mb(tmp1(1)),sum)
8675         call Pack_cc_dot(1,dcpl_mb(tmp1(1)),dcpl_mb(tmp1(1)),sum)
8676
8677         error = error + sum
8678      end do
8679      call D1dB_SumAll(error)
8680      error = dsqrt(error)/size
8681
8682      value = BA_pop_stack(tmp1(2))
8683      if (.not. value)
8684     > call errquit('error popping stack memory in psi_error',0, MA_ERR)
8685
8686
8687      psi_error = error
8688      return
8689      end
8690
8691*     ***************************
8692*     *				*
8693*     *		rho_error	*
8694*     *				*
8695*     ***************************
8696      real*8 function rho_error()
8697      implicit none
8698
8699#include "errquit.fh"
8700#include "bafdecls.fh"
8701#include "psi.fh"
8702
8703*     ***** local variables ****
8704      logical value
8705      integer k,nx,ny,nz
8706      real*8  error,scale
8707      integer tmp1(2)
8708
8709      real*8 e1,e2
8710      common /eenergy_tmp_common/ e1,e2
8711
8712*     ***** external functions *****
8713      real*8   lattice_omega
8714      external lattice_omega
8715
8716      value = BA_push_get(mt_dbl,(2*nfft3d),'tmp1',tmp1(2),tmp1(1))
8717      if (.not. value)
8718     >   call errquit('out of stack memory in rho_error',0, MA_ERR)
8719
8720
8721      call D3dB_nx(1,nx)
8722      call D3dB_ny(1,ny)
8723      call D3dB_nz(1,nz)
8724      scale = lattice_omega()
8725
8726      scale = (scale)/dble(nx*ny*nz)
8727*     scale = (scale)/dble(nx*ny*nz)
8728*     scale = (scale*scale)
8729
8730!$OMP DO
8731      do k=1,(2*nfft3d)
8732         dbl_mb(tmp1(1)+k-1) = (dbl_mb(rho2(1)+k-1)
8733     >                         -dbl_mb(rho1(1)+k-1))
8734         dbl_mb(tmp1(1)+k-1) = dbl_mb(tmp1(1)+k-1)
8735     >                      + (dbl_mb(rho2(1)+k-1+(ispin-1)*(2*nfft3d))
8736     >                        -dbl_mb(rho1(1)+k-1+(ispin-1)*(2*nfft3d)))
8737      end do
8738!$OMP END DO
8739      call D3dB_rr_dot(1,dbl_mb(tmp1(1)),dbl_mb(tmp1(1)),e1)
8740      error = e1*scale
8741*     error = dsqrt(error)
8742
8743
8744      value = BA_pop_stack(tmp1(2))
8745      if (.not. value)
8746     > call errquit('error popping stack memory in rho_error',0, MA_ERR)
8747
8748
8749      rho_error = error
8750      return
8751      end
8752
8753
8754*     ***************************
8755*     *                         *
8756*     *         psi_a_sum       *
8757*     *                         *
8758*     ***************************
8759      real*8 function psi_a_sum(npack1,psi)
8760      implicit none
8761      integer npack1
8762      complex*16 psi(*)
8763
8764      integer k
8765      real*8 a,tmp
8766
8767      a = 0.0d0
8768      do k=1,npack1
8769         tmp = dble(psi(k))
8770         a = a + tmp*tmp
8771      end do
8772      call D3dB_SumAll(a)
8773
8774      psi_a_sum = a
8775      return
8776      end
8777
8778
8779
8780
8781*     ***************************
8782*     *                         *
8783*     *         psi_b_sum       *
8784*     *                         *
8785*     ***************************
8786      real*8 function psi_b_sum(npack1,psi)
8787      implicit none
8788      integer npack1
8789      complex*16 psi(*)
8790
8791
8792      integer k
8793      real*8 b,tmp
8794
8795      b = 0.0d0
8796      do k=1,npack1
8797         tmp = dimag(psi(k))
8798         b = b + tmp*tmp
8799      end do
8800      call D3dB_SumAll(b)
8801
8802      psi_b_sum = b
8803      return
8804      end
8805
8806*     **************************************
8807*     *                                    *
8808*     *          psi_symm_project          *
8809*     *                                    *
8810*     **************************************
8811      subroutine psi_a_project(npack1,psi)
8812      implicit none
8813      integer npack1
8814      complex*16 psi(*)
8815      integer k
8816      real*8 tmp
8817      do k=1,npack1
8818        tmp    = dble(psi(k))
8819        psi(k) = dcmplx(tmp,0.0d0)
8820      end do
8821      return
8822      end
8823      subroutine psi_b_project(npack1,psi)
8824      implicit none
8825      integer npack1
8826      complex*16 psi(*)
8827      integer k
8828      real*8 tmp
8829      do k=1,npack1
8830        tmp    = dimag(psi(k))
8831        psi(k) = dcmplx(0.0d0,tmp)
8832      end do
8833      return
8834      end
8835
8836      subroutine psi_symm_project(ispin,ne,npack1,psi1)
8837      implicit none
8838      integer ispin,ne(2),npack1
8839      complex*16 psi1(npack1,*)
8840
8841      integer i
8842      real*8   a,b
8843      real*8   psi_a_sum,psi_b_sum
8844      external psi_a_sum,psi_b_sum
8845
8846      do i=1,(ne(1)+ne(2))
8847          a = psi_a_sum(npack1,psi1(1,i))
8848          b = psi_b_sum(npack1,psi1(1,i))
8849          if (a.ge.b) then
8850             call psi_a_project(npack1,psi1(1,i))
8851          else
8852             call psi_b_project(npack1,psi1(1,i))
8853          end if
8854      end do
8855      return
8856      end
8857
8858      subroutine psi_ab_gen_irrep_names(virtual)
8859      implicit none
8860      logical virtual
8861
8862#include "bafdecls.fh"
8863#include "errquit.fh"
8864#include "psi.fh"
8865
8866      integer irreps(2)
8867      common / ab_irrep / irreps
8868
8869      integer k,n,ptr,nn
8870      real*8  a,b,tmpa,tmpb
8871
8872      if (virtual) then
8873          ptr = psi1_excited(1)
8874          nn  = ne_excited(1)+ne_excited(2)
8875      else
8876         ptr = psi1(1)
8877         nn  = ne(1)+ne(2)
8878      end if
8879
8880      if (.not.BA_alloc_get(mt_int,nn,
8881     >                     'irreps',irreps(2),irreps(1)))
8882     > call errquit('psi_ab_gen_irrep_names',0, MA_ERR)
8883
8884      do n=1,nn
8885         a = 0.0d0
8886         b = 0.0d0
8887         do k=1,npack1
8888            tmpa = dble( dcpl_mb(ptr+k-1+(n-1)*npack1))
8889            tmpb = dimag(dcpl_mb(ptr+k-1+(n-1)*npack1))
8890            a = a + tmpa*tmpa
8891            b = b + tmpb*tmpb
8892         end do
8893         call D3dB_SumAll(a)
8894         call D3dB_SumAll(b)
8895
8896         if      ((b .lt. 1.0d-6).and.(a .gt. 1.0d-6)) then
8897            int_mb(irreps(1)+n-1) = 1
8898         else if ((a .lt. 1.0d-6).and.(b .gt. 1.0d-6)) then
8899            int_mb(irreps(1)+n-1) = -1
8900         else
8901            int_mb(irreps(1)+n-1) = 0
8902         end if
8903      end do
8904
8905
8906      return
8907      end
8908
8909      subroutine psi_ab_kill_irrep_names()
8910      implicit none
8911
8912#include "bafdecls.fh"
8913#include "errquit.fh"
8914
8915      integer irreps(2)
8916      common / ab_irrep / irreps
8917
8918      if (.not.BA_free_heap(irreps(2)))
8919     >  call errquit('psi_ab_gen_irrep_names: error freeing heap',
8920     >               0, MA_ERR)
8921
8922      return
8923      end
8924
8925
8926
8927*     **************************************
8928*     *                                    *
8929*     *         psi_ab_irrep_name          *
8930*     *                                    *
8931*     **************************************
8932
8933*     This function resturns
8934*        '[ag]' - if psi(n) is purely real
8935*        '[au]' - if psi(n) is purely imaginary
8936*        '    ' - if psi(n) is mixed
8937*
8938*   Not psi_ab_gen_irrep_names needs to be called before this is used.
8939*
8940      character*4 function psi_ab_irrep_name(n)
8941      implicit none
8942      integer n
8943
8944#include "bafdecls.fh"
8945
8946      integer irreps(2)
8947      common / ab_irrep / irreps
8948
8949      character*4 abvalue
8950
8951      if      (int_mb(irreps(1)+n-1).eq.1) then
8952         abvalue = '[ag]'
8953      else if (int_mb(irreps(1)+n-1).eq.-1) then
8954         abvalue = '[au]'
8955      else
8956         abvalue = '    '
8957      end if
8958
8959      psi_ab_irrep_name = abvalue
8960      return
8961      end
8962
8963
8964*     ***************************
8965*     *                         *
8966*     *   psi1_crystal_dipole   *
8967*     *                         *
8968*     ***************************
8969*
8970*     Uses - electron_crystal_dipole
8971*
8972      subroutine psi1_crystal_dipole(dipole)
8973      implicit none
8974      real*8 dipole(3)
8975
8976#include "bafdecls.fh"
8977#include "psi.fh"
8978
8979      call Calculate_Resta_Dipole(.true.,ispin,ne,neq,npack1,nfft3d,
8980     >                            dcpl_mb(psi1(1)),dipole)
8981
8982      return
8983      end
8984
8985
8986
8987*     ***************************
8988*     *				*
8989*     *		rho_dipole	*
8990*     *				*
8991*     ***************************
8992*
8993*     Uses - Calculate_dipole (pspw/lib/psi/dipole.f)
8994*
8995      subroutine rho_dipole(dipole)
8996      implicit none
8997      real*8 dipole(3)
8998
8999#include "bafdecls.fh"
9000#include "psi.fh"
9001
9002      call Calculate_Dipole(ispin,ne,2*nfft3d,dbl_mb(rho1(1)),dipole)
9003      return
9004      end
9005
9006
9007*     ***************************
9008*     *				*
9009*     *		psi_ispin	*
9010*     *				*
9011*     ***************************
9012      integer function psi_ispin()
9013      implicit none
9014
9015#include "bafdecls.fh"
9016#include "psi.fh"
9017
9018      psi_ispin = ispin
9019      return
9020      end
9021
9022
9023*     ***************************
9024*     *				*
9025*     *		psi_ne		*
9026*     *				*
9027*     ***************************
9028      integer function psi_ne(ms)
9029      implicit none
9030      integer ms
9031
9032#include "bafdecls.fh"
9033#include "psi.fh"
9034
9035      psi_ne = ne(ms)
9036      return
9037      end
9038
9039*     ***************************
9040*     *				*
9041*     *		psi_neq		*
9042*     *				*
9043*     ***************************
9044      integer function psi_neq(ms)
9045      implicit none
9046      integer ms
9047
9048#include "bafdecls.fh"
9049#include "psi.fh"
9050
9051      psi_neq = neq(ms)
9052      return
9053      end
9054
9055
9056
9057*     ***************************
9058*     *                         *
9059*     *    psi_cpmd_start       *
9060*     *                         *
9061*     ***************************
9062      subroutine psi_cpmd_start()
9063      implicit none
9064
9065#include "bafdecls.fh"
9066#include "errquit.fh"
9067#include "psi.fh"
9068
9069      logical value
9070
9071      value = BA_alloc_get(mt_dbl,2*ispin*nfft3d,'rho0',rho0(2),rho0(1))
9072      if (.not.value)
9073     >   call errquit('psi_cpmd_start',0,MA_ERR)
9074
9075      call Parallel_shared_vector_copy(.true.,2*ispin*nfft3d,
9076     >                                 dbl_mb(rho1(1)),dbl_mb(rho0(1)))
9077      return
9078      end
9079
9080*     ***************************
9081*     *                         *
9082*     *    psi_cpmd_end         *
9083*     *                         *
9084*     ***************************
9085      subroutine psi_cpmd_end()
9086      implicit none
9087
9088#include "bafdecls.fh"
9089#include "errquit.fh"
9090#include "psi.fh"
9091
9092      if (.not.BA_free_heap(rho0(2)))
9093     >   call errquit('psi_cpmd_end',0,MA_ERR)
9094      return
9095      end
9096
9097
9098*     ***************************
9099*     *                         *
9100*     *    psi_cpmd_step        *
9101*     *                         *
9102*     ***************************
9103      subroutine psi_cpmd_step(dte)
9104      implicit none
9105      real*8 dte
9106
9107#include "bafdecls.fh"
9108#include "errquit.fh"
9109#include "psi.fh"
9110
9111      logical  control_precondition
9112      external control_precondition
9113      integer  control_ks_algorithm
9114      external control_ks_algorithm
9115      real*8   control_tole
9116      external control_tole
9117
9118
9119*     **** psi2 = 2*psi1 - psi0 + dt*dt/fmass*Hpsi ****
9120c      call electron_cpmd_update(dcpl_mb(psi0(1)),
9121c     >                          dcpl_mb(psi1(1)),
9122c     >                          dcpl_mb(psi2(1)),
9123c     >                          dbl_mb(hml(1)),
9124c     >                          dte)
9125c      call Dneall_f_ortho(0,dcpl_mb(psi2(1)),npack1)
9126c      write(*,*) "psi1 ortho:"
9127c      call OrthoCheck_geo(ispin,ne,dcpl_mb(psi1(1)))
9128c*     **** lagrange multiplier corrections ****
9129c      if (pspw_SIC().or.occupation_on) then
9130c        call psi_lmbda_sic(ispin,ne,(neq(1)+neq(2)),npack1,
9131c     >                 dcpl_mb(psi1(1)),dcpl_mb(psi2(1)),dte,
9132c     >                 dbl_mb(lmd_cpmd(1)),
9133c     >                 dbl_mb(tmp_L_cpmd(1)),ierr)
9134c      else
9135c        call psi_lmbda(ispin,ne,(neq(1)+neq(2)),npack1,
9136c     >                 dcpl_mb(psi1(1)),dcpl_mb(psi2(1)),dte,
9137c     >                 dbl_mb(lmd_cpmd(1)),
9138c     >                 dbl_mb(tmp_L_cpmd(1)),ierr)
9139c      end if
9140c      write(*,*) "psi2 ortho:"
9141c      call OrthoCheck_geo(ispin,ne,dcpl_mb(psi2(1)))
9142
9143      call Parallel_shared_vector_copy(.true.,2*ispin*nfft3d,
9144     >                          dbl_mb(rho1(1)),dbl_mb(rho2(1)))
9145      call DSCAL_OMP(2*ispin*nfft3d, 2.0d0,dbl_mb(rho2(1)),1)
9146      call DAXPY_OMP(2*ispin*nfft3d,-1.0d0,
9147     >               dbl_mb(rho0(1)),1,dbl_mb(rho2(1)),1)
9148      call Parallel_shared_vector_copy(.true.,2*ispin*nfft3d,
9149     >                          dbl_mb(rho1(1)),dbl_mb(rho0(1)))
9150
9151      call psi_set_density(1,dbl_mb(rho2(1)))
9152
9153*     **** diaganolize KS matrix ****
9154      call psi_KS_update(1,
9155     >                   control_ks_algorithm(),
9156     >                   control_precondition(),
9157     >                   control_tole())
9158
9159      return
9160      end
9161
9162
9163
9164*     ***************************
9165*     *				*
9166*     *	    psi_initialize 	*
9167*     *				*
9168*     ***************************
9169
9170      logical function psi_initialize()
9171      implicit none
9172
9173
9174#include "bafdecls.fh"
9175#include "btdb.fh"
9176#include "stdio.fh"
9177#include "util.fh"
9178#include "errquit.fh"
9179#include "psi.fh"
9180
9181*     ***** rhoall common block ****
9182      integer rho1_all(2)
9183      integer rho2_all(2)
9184      common / rhoall_block / rho1_all,rho2_all
9185
9186*     **** local variables ****
9187      integer MASTER,taskid
9188      parameter (MASTER=0)
9189      logical value,psi_nogrid
9190      integer nemax
9191      real*8 sum1,sum2,sum3
9192      integer hversion,hnfft(3),hispin,hne(2)
9193      real*8 hunita(3,3)
9194      integer rtdb,ind,vers
9195      integer  control_rtdb,control_ngrid,control_symmetry
9196      external control_rtdb,control_ngrid,control_symmetry
9197      integer  psi_get_version
9198      external psi_get_version
9199      character*50 filename
9200      character*50 control_input_psi
9201      external     control_input_psi
9202      logical  wvfnc_expander,Dneall_m_allocate,band_reformat_c_wvfnc
9203      external wvfnc_expander,Dneall_m_allocate,band_reformat_c_wvfnc
9204      logical  psp_pawexist,control_print
9205      external psp_pawexist,control_print
9206      integer          control_fractional_smeartype
9207      double precision control_fractional_kT
9208      external         control_fractional_smeartype
9209      external         control_fractional_kT
9210      logical          control_ortho
9211      external         control_ortho
9212
9213
9214      ne_excited(1) = 0
9215      ne_excited(2) = 0
9216
9217*     **** reformat wavefunction if it is a band wavefunction ****
9218      vers = psi_get_version()
9219      if (vers.eq.5) then
9220        call Parallel_taskid(taskid)
9221        if (taskid.eq.MASTER) then
9222          value= band_reformat_c_wvfnc(1)
9223        end if
9224      end if
9225      pawexist = psp_pawexist()
9226
9227
9228*     *****  get ispin,ne,neq,nfft3d,npack0,npack1 ****
9229      call psi_get_ne_occupation(ispin,ne,smearoccupation)
9230      call Dneall_neq(neq)
9231      call D3dB_nfft3d(1,nfft3d)
9232      call Pack_npack(1,npack1)
9233      call Pack_npack(0,npack0)
9234      nemax = ne(1)+ne(2)
9235      occupation_on = .false.
9236      if (smearoccupation.gt.0) occupation_on = .true.
9237
9238
9239*     **** allocate memory ****
9240      value = BA_alloc_get(mt_dcpl,npack1*(neq(1)+neq(2)),
9241     >                     'psi2',psi2(2),psi2(1))
9242      value = value.and.
9243     >        BA_alloc_get(mt_dcpl,npack1*(neq(1)+neq(2)),
9244     >                     'psi1',psi1(2),psi1(1))
9245      if (pawexist) then
9246      value = value.and.
9247     >        BA_alloc_get(mt_dcpl,npack1*(neq(1)+neq(2)),
9248     >                     'spsi1',spsi1(2),spsi1(1))
9249      end if
9250      value = value.and.
9251     >        BA_alloc_get(mt_dbl,4*nfft3d,
9252     >                     'rho1',rho1(2),rho1(1))
9253      value = value.and.
9254     >        BA_alloc_get(mt_dbl,4*nfft3d,
9255     >                     'rho2',rho2(2),rho2(1))
9256      value = value.and.
9257     >        BA_alloc_get(mt_dcpl,npack0,
9258     >                     'dng1',dng1(2),dng1(1))
9259      value = value.and.
9260     >        BA_alloc_get(mt_dcpl,npack0,
9261     >                     'dng2',dng2(2),dng2(1))
9262c      value = value.and.
9263c     >        BA_alloc_get(mt_dbl,(2*nemax*nemax),'hml',hml(2),hml(1))
9264      value = value.and.Dneall_m_allocate(0,hml)
9265
9266      value = value.and.
9267     >        BA_alloc_get(mt_dbl,(2*nemax),'eig',eig(2),eig(1))
9268
9269      if (occupation_on) then
9270        value = value.and.
9271     >        BA_alloc_get(mt_dbl,(2*nemax),'occ1',occ1(2),occ1(1))
9272        value = value.and.
9273     >        BA_alloc_get(mt_dbl,(2*nemax),'occ2',occ2(2),occ2(1))
9274        smeartype = control_fractional_smeartype()
9275        smearkT   = control_fractional_kT()
9276      end if
9277
9278      value = value.and.
9279     >        BA_alloc_get(mt_dbl,4*nfft3d,
9280     >                     'rho1_all',rho1_all(2),rho1_all(1))
9281      value = value.and.
9282     >        BA_alloc_get(mt_dbl,4*nfft3d,
9283     >                     'rho2_all',rho2_all(2),rho2_all(1))
9284      if (.not. value) call errquit('out of heap memory',0, MA_ERR)
9285c      call dcopy(4*nfft3d,0.0d0,0,dbl_mb(rho1_all(1)),1)
9286c      call dcopy(4*nfft3d,0.0d0,0,dbl_mb(rho2_all(1)),1)
9287      call Parallel_shared_vector_zero(.true.,4*nfft3d,
9288     >                    dbl_mb(rho1_all(1)))
9289      call Parallel_shared_vector_zero(.true.,4*nfft3d,
9290     >                    dbl_mb(rho2_all(1)))
9291
9292*     *****  read initial wavefunctions into psi1  ****
9293      rtdb = control_rtdb()
9294      if (.not.btdb_get(rtdb,'nwpw:psi_nogrid',
9295     >                  mt_log,1,psi_nogrid))
9296     >   psi_nogrid = .true.
9297
9298      if (psi_nogrid) then
9299
9300        call psi_get_header(hversion,hnfft,hunita,hispin,hne)
9301
9302        if ( (hnfft(1).ne.control_ngrid(1)) .or.
9303     >       (hnfft(2).ne.control_ngrid(2)) .or.
9304     >       (hnfft(3).ne.control_ngrid(3)) ) then
9305
9306        hnfft(1) = control_ngrid(1)
9307        hnfft(2) = control_ngrid(2)
9308        hnfft(3) = control_ngrid(3)
9309        call Parallel_taskid(taskid)
9310
9311        call ga_sync()
9312        value = btdb_parallel(.false.)
9313        call ga_sync()
9314        if (taskid.eq.MASTER) then
9315
9316          filename =  control_input_psi()
9317
9318          ind = index(filename,' ') - 1
9319          if (.not. btdb_cput(rtdb,'xpndr:old_wavefunction_filename',
9320     >                    1,filename(1:ind)))
9321     >     call errquit(
9322     >     'wvfnc_expander_input: btdb_cput failed', 0, RTDB_ERR)
9323
9324          if (.not. btdb_cput(rtdb,'xpndr:new_wavefunction_filename',
9325     >                    1,filename(1:ind)))
9326     >     call errquit(
9327     >     'wvfnc_expander_input: btdb_cput failed', 0, RTDB_ERR)
9328
9329          if (.not. btdb_put(rtdb,'xpndr:ngrid',mt_int,3,hnfft))
9330     >     call errquit(
9331     >     'wvfnc_expander_input: btdb_put failed', 0, RTDB_ERR)
9332
9333          if (control_print(print_medium)) then
9334            write(luout,*)
9335            write(luout,*) "Grid is being converted:"
9336            write(luout,*) "------------------------"
9337            write(luout,*)
9338            write(luout,*) "To turn off automatic grid conversion:"
9339            write(luout,*)
9340            write(luout,*) "set nwpw:psi_nogrid .false."
9341            write(luout,*)
9342          endif
9343          value = wvfnc_expander(rtdb)
9344
9345        end if
9346        call ga_sync()
9347        value = btdb_parallel(.true.)
9348        value = .true.
9349
9350      end if
9351
9352      end if
9353
9354      call psi_read(ispin,ne,dcpl_mb(psi1(1)),
9355     >              smearoccupation,dbl_mb(occ1(1)))
9356
9357      if (occupation_on) then
9358         call frac_occ_set(rtdb,ispin,ne,dbl_mb(occ1(1)))
9359      end if
9360
9361      call psi_history_read(ispin,ne,
9362     >                      dcpl_mb(psi1(1)),
9363     >                      dcpl_mb(psi2(1)))
9364
9365
9366*     **** force inversion symmetry ****
9367      if (control_symmetry().eq.1)  then
9368         call Parallel_taskid(taskid)
9369         if ((taskid.eq.MASTER).and.
9370     >       (control_print(print_medium))) then
9371         write(luout,*)
9372         write(luout,*)
9373     >   "Projecting wavefunctions to have inversion symmetry"
9374         write(luout,*)
9375         end if
9376         call psi_symm_project(ispin,neq,npack1,dcpl_mb(psi1(1)))
9377      end if
9378
9379*     **** Ortho Check ****
9380      call psi_1ortho_check_fix()
9381c      if (pawexist) then
9382c         call psp_overlap_S(ispin,neq,dcpl_mb(psi1(1)),dcpl_mb(psi2(1)))
9383c         call Grsm_gg_trace(npack1,(neq(1)+neq(2)),
9384c     >                        dcpl_mb(psi1(1)),
9385c     >                        dcpl_mb(psi2(1)),
9386c     >                        sum2)
9387c      else
9388c         call Grsm_gg_trace(npack1,(neq(1)+neq(2)),
9389c     >                        dcpl_mb(psi1(1)),
9390c     >                        dcpl_mb(psi1(1)),
9391c     >                        sum2)
9392c      end if
9393c      call D1dB_SumAll(sum2)
9394c
9395c
9396c      sum1 = dble(ne(1) + ne(2))
9397c      if ((control_ortho()).and.(dabs(sum2-sum1).gt.1.0d-10)) then
9398c         call Parallel_taskid(taskid)
9399c         if (pawexist) then
9400c         call Dneall_f_Sortho(0,dcpl_mb(psi1(1)),
9401c     >                          dcpl_mb(psi2(1)),npack1)
9402c         call psp_overlap_S(ispin,neq,dcpl_mb(psi1(1)),
9403c     >                                dcpl_mb(psi2(1)))
9404c         call Grsm_gg_trace(npack1,(neq(1)+neq(2)),
9405c     >                        dcpl_mb(psi1(1)),
9406c     >                        dcpl_mb(psi2(1)),
9407c     >                        sum3)
9408c         else
9409cc         call Dneall_f_ortho(0,dcpl_mb(psi1(1)),npack1)
9410c         call Dneall_f_GramSchmidt(0,dcpl_mb(psi1(1)),npack1)
9411c         call Grsm_gg_trace(npack1,(neq(1)+neq(2)),
9412c     >                        dcpl_mb(psi1(1)),
9413c     >                        dcpl_mb(psi1(1)),
9414c     >                        sum3)
9415c         end if
9416c         call D1dB_SumAll(sum3)
9417c         if ((taskid.eq.MASTER).and.(control_print(print_medium)))
9418c     >    write(luout,*)
9419c     >     "Warning - Gram-Schmidt being performed on psi:",
9420c     >               sum1,sum2,sum3,dabs(sum2-sum1)
9421c
9422c
9423c      end if
9424
9425      psi_initialize = value
9426      return
9427      end
9428
9429*     ***************************
9430*     *				*
9431*     *	 psi_1ortho_check_fix   *
9432*     *				*
9433*     ***************************
9434      subroutine psi_1ortho_check_fix()
9435      implicit none
9436
9437#include "bafdecls.fh"
9438#include "stdio.fh"
9439#include "util.fh"
9440#include "psi.fh"
9441
9442*     **** local variables ****
9443      integer taskid,MASTER
9444      parameter (MASTER=0)
9445
9446      real*8 sum2,sum1,sum3
9447
9448*     **** external functions ****
9449      logical  control_ortho,control_print
9450      external control_ortho,control_print
9451
9452*     **** Ortho Check ****
9453      if (pawexist) then
9454         call psp_overlap_S(ispin,neq,dcpl_mb(psi1(1)),dcpl_mb(psi2(1)))
9455         call Grsm_gg_trace(npack1,(neq(1)+neq(2)),
9456     >                        dcpl_mb(psi1(1)),
9457     >                        dcpl_mb(psi2(1)),
9458     >                        sum2)
9459      else
9460         call Grsm_gg_trace(npack1,(neq(1)+neq(2)),
9461     >                        dcpl_mb(psi1(1)),
9462     >                        dcpl_mb(psi1(1)),
9463     >                        sum2)
9464      end if
9465      call D1dB_SumAll(sum2)
9466
9467
9468      sum1 = dble(ne(1) + ne(2))
9469      if ((control_ortho()).and.(dabs(sum2-sum1).gt.1.0d-10)) then
9470         call Parallel_taskid(taskid)
9471         if (pawexist) then
9472         call Dneall_f_Sortho(0,dcpl_mb(psi1(1)),
9473     >                          dcpl_mb(psi2(1)),npack1)
9474         call psp_overlap_S(ispin,neq,dcpl_mb(psi1(1)),
9475     >                                dcpl_mb(psi2(1)))
9476         call Grsm_gg_trace(npack1,(neq(1)+neq(2)),
9477     >                        dcpl_mb(psi1(1)),
9478     >                        dcpl_mb(psi2(1)),
9479     >                        sum3)
9480         else
9481c         call Dneall_f_ortho(0,dcpl_mb(psi1(1)),npack1)
9482         call Dneall_f_GramSchmidt(0,dcpl_mb(psi1(1)),npack1)
9483         call Grsm_gg_trace(npack1,(neq(1)+neq(2)),
9484     >                        dcpl_mb(psi1(1)),
9485     >                        dcpl_mb(psi1(1)),
9486     >                        sum3)
9487         end if
9488         call D1dB_SumAll(sum3)
9489         if ((taskid.eq.MASTER).and.(control_print(print_medium))) then
9490             write(luout,321) sum1,sum2,sum3,dabs(sum2-sum1)
9491         end if
9492 321     format(/" Warning - K.S. orbitals are not orthonormal. ",
9493     >      "Applying Gram-Schmidt orthonormalization."/,
9494     >       "         - exact norm=",E12.6," norm=",E12.6,
9495     >       " corrected norm=",E12.6,
9496     >       " (error=",E12.6,")"/)
9497
9498      end if
9499      return
9500      end
9501
9502
9503
9504*     ***************************
9505*     *				*
9506*     *	  psi_tmp_write  	*
9507*     *				*
9508*     ***************************
9509      subroutine psi_tmp_write()
9510      implicit none
9511
9512#include "bafdecls.fh"
9513#include "psi.fh"
9514
9515*     ***** write psi1 wavefunctions ****
9516      call psi_write(ispin,ne,dcpl_mb(psi1(1)),
9517     >               smearoccupation,dbl_mb(occ1(1)))
9518
9519      return
9520      end
9521
9522
9523*     ************************************
9524*     *				         *
9525*     *	    psi_tmp_write_full_filename  *
9526*     *				         *
9527*     ************************************
9528      subroutine psi_tmp_write_full_filename(full_filename)
9529      implicit none
9530      character*(*) full_filename
9531
9532#include "bafdecls.fh"
9533#include "psi.fh"
9534
9535*     ***** write psi1 wavefunctions ****
9536      call psi_write_full_filename(full_filename,
9537     >                             ispin,ne,dcpl_mb(psi1(1)),
9538     >                             smearoccupation,dbl_mb(occ1(1)))
9539
9540      return
9541      end
9542
9543*     ************************************
9544*     *				         *
9545*     *	    psi_tmp_read_full_filename   *
9546*     *				         *
9547*     ************************************
9548      subroutine psi_tmp_read_full_filename(full_filename)
9549      implicit none
9550      character*(*) full_filename
9551
9552#include "bafdecls.fh"
9553#include "psi.fh"
9554
9555      call psi_read_full_filename(full_filename,
9556     >                            ispin,ne,dcpl_mb(psi1(1)),
9557     >                            smearoccupation,dbl_mb(occ1(1)))
9558      return
9559      end
9560
9561
9562
9563*     ***************************
9564*     *				*
9565*     *		psi_finalize	*
9566*     *				*
9567*     ***************************
9568
9569      logical function psi_finalize(wpsi)
9570      implicit none
9571      logical wpsi
9572
9573#include "errquit.fh"
9574#include "bafdecls.fh"
9575#include "psi.fh"
9576
9577*     ***** rhoall common block ****
9578      integer rho1_all(2)
9579      integer rho2_all(2)
9580      common / rhoall_block / rho1_all,rho2_all
9581
9582
9583*     **** local variables ****
9584      logical value
9585
9586      logical  Dneall_m_free
9587      external Dneall_m_free
9588
9589*     ***** write psi1 wavefunctions ****
9590      if (wpsi) then
9591        call psi_write(ispin,ne,dcpl_mb(psi1(1)),
9592     >                 smearoccupation,dbl_mb(occ1(1)))
9593        call psi_history_write(ispin,ne,dcpl_mb(psi1(1)))
9594      end if
9595
9596      value = BA_free_heap(eig(2))
9597      value = value.and.Dneall_m_free(hml)
9598      value = value.and.BA_free_heap(dng2(2))
9599      value = value.and.BA_free_heap(dng1(2))
9600      value = value.and.BA_free_heap(rho2(2))
9601      value = value.and.BA_free_heap(rho1(2))
9602      value = value.and.BA_free_heap(psi2(2))
9603      value = value.and.BA_free_heap(psi1(2))
9604      if (pawexist) then
9605          value = value.and.BA_free_heap(spsi1(2))
9606      end if
9607      value = value.and.BA_free_heap(rho2_all(2))
9608      value = value.and.BA_free_heap(rho1_all(2))
9609      if (occupation_on) then
9610         value = value.and.BA_free_heap(occ2(2))
9611         value = value.and.BA_free_heap(occ1(2))
9612      end if
9613
9614      if (.not. value)
9615     >  call errquit('psi_finalize: error freeing heap',0, MA_ERR)
9616
9617      psi_finalize = value
9618      return
9619      end
9620
9621
9622*     ***************************
9623*     *                         *
9624*     *      psi_ne_excited     *
9625*     *                         *
9626*     ***************************
9627      integer function psi_ne_excited(ms)
9628      implicit none
9629      integer ms
9630
9631#include "bafdecls.fh"
9632#include "psi.fh"
9633
9634      psi_ne_excited = ne_excited(ms)
9635      return
9636      end
9637
9638*     ***************************
9639*     *				*
9640*     *   psi_2epsi_gradients	*
9641*     *				*
9642*     ***************************
9643      logical function psi_2epsi_gradients(gn)
9644      implicit none
9645      integer gn
9646
9647#include "bafdecls.fh"
9648#include "btdb.fh"
9649#include "psi.fh"
9650#include "errquit.fh"
9651
9652*     **** local variables ****
9653      logical value
9654      integer nfac,nex(2),nemax,nexmax,ii,n
9655      real*8 tsum,tsum0
9656
9657      nfac = 3
9658      if (gn.eq.2) then
9659         nfac = 9
9660      else if (gn.eq.3) then
9661         nfac = 19
9662      else
9663         nfac = gn
9664      end if
9665      nex(1) = ne(1)*nfac
9666      nex(2) = ne(2)*nfac
9667      nemax   = ne(1) + ne(2)
9668      nexmax  = nex(1) + nex(2)
9669
9670*     **** allocate memory ****
9671      value = BA_alloc_get(mt_dcpl,npack1*(nexmax),
9672     >         'psi1_excited',psi1_excited(2),psi1_excited(1))
9673      if (.not.value)
9674     >  call errquit('psi_2epsi_gradients: out of heap memory',0,MA_ERR)
9675
9676      do ii=1,nfac
9677         do n=1,nemax
9678            call Pack_cc_multiplegradients(1,ii,
9679     >             dcpl_mb(psi1(1)+(n-1)*npack1),
9680     >             dcpl_mb(psi1_excited(1)+((n-1)+(ii-1)*nemax)*npack1))
9681            call Pack_cc_dot(1,
9682     >             dcpl_mb(psi1_excited(1)+((n-1)+(ii-1)*nemax)*npack1),
9683     >             dcpl_mb(psi1_excited(1)+((n-1)+(ii-1)*nemax)*npack1),
9684     >             tsum)
9685            call Pack_cc_dot(1,
9686     >             dcpl_mb(psi1(1)+(n-1)*npack1),
9687     >             dcpl_mb(psi1(1)+(n-1)*npack1),
9688     >             tsum0)
9689            write(*,*) "ii,n,tsum=",ii,n,tsum,tsum0
9690         end do
9691      end do
9692
9693      call epsi_write(ispin,nex,dcpl_mb(psi1_excited(1)))
9694      value = BA_free_heap(psi1_excited(2))
9695      if (.not.value)
9696     >  call errquit('psi_2epsi_gradients:freeing heap',0,MA_ERR)
9697
9698      return
9699      end
9700
9701
9702*     ***************************
9703*     *				*
9704*     *     epsi_initialize 	*
9705*     *				*
9706*     ***************************
9707
9708      logical function epsi_initialize()
9709      implicit none
9710
9711#include "bafdecls.fh"
9712#include "btdb.fh"
9713#include "psi.fh"
9714#include "errquit.fh"
9715
9716*     **** local variables ****
9717      integer MASTER,taskid
9718      parameter (MASTER=0)
9719      logical value,psi_nogrid
9720      integer nemax,ispin0
9721
9722      integer hversion,hnfft(3),hispin,hne(2)
9723      real*8 hunita(3,3),sum1,sum2
9724      integer rtdb,ind,nn
9725      integer  control_rtdb,control_ngrid
9726      external control_rtdb,control_ngrid
9727      character*50 filename
9728      character*50 control_input_epsi
9729      external     control_input_epsi
9730      logical  wvfnc_expander
9731      external wvfnc_expander
9732      integer  control_symmetry
9733      external control_symmetry
9734
9735
9736*     ***** get ispin, and ne, and nfft3d ****
9737      call psi_get_ne_excited(ispin0,ne_excited)
9738      nemax  = ne_excited(1)  + ne_excited(2)
9739      nn = ne_excited(1)*ne_excited(1) + ne_excited(2)*ne_excited(2)
9740
9741*     **** allocate memory ****
9742      value = BA_alloc_get(mt_dcpl,npack1*(nemax),
9743     >         'psi2_excited',psi2_excited(2),psi2_excited(1))
9744      value = value.and.
9745     >        BA_alloc_get(mt_dcpl,npack1*(nemax),
9746     >         'psi1_excited',psi1_excited(2),psi1_excited(1))
9747      value = value.and.
9748     >        BA_alloc_get(mt_dbl,(2*nemax),'eig_excited',
9749     >                     eig_excited(2),eig_excited(1))
9750      value = value.and.
9751     >        BA_alloc_get(mt_dbl,nn,'hml_excited',
9752     >                     hml_excited(2),hml_excited(1))
9753      if (.not.value)
9754     >   call errquit('epsi_initialize: out of heap memory',0,MA_ERR)
9755
9756      !call dcopy(2*npack1*nemax,0.0d0,0,dcpl_mb(psi2_excited(1)),1)
9757      !call dcopy(2*npack1*nemax,0.0d0,0,dcpl_mb(psi1_excited(1)),1)
9758      !call dcopy(2*nemax,0.0d0,0,dbl_mb(eig_excited(1)),1)
9759      !call dcopy(nn,0.0d0,0,dbl_mb(hml_excited(1)),1)
9760      call Parallel_shared_vector_zero(.true.,2*npack1*nemax,
9761     >                    dcpl_mb(psi2_excited(1)))
9762      call Parallel_shared_vector_zero(.true.,2*npack1*nemax,
9763     >                    dcpl_mb(psi1_excited(1)))
9764      call Parallel_shared_vector_zero(.true.,2*nemax,
9765     >                    dbl_mb(eig_excited(1)))
9766      call Parallel_shared_vector_zero(.true.,nn,
9767     >                    dbl_mb(hml_excited(1)))
9768
9769
9770*     *****  read initial wavefunctions into psi1  ****
9771      rtdb = control_rtdb()
9772      if (.not.btdb_get(rtdb,'nwpw:psi_nogrid',
9773     >                  mt_log,1,psi_nogrid))
9774     >   psi_nogrid = .true.
9775
9776      if (psi_nogrid) then
9777
9778
9779        filename =  control_input_epsi()
9780        call psi_get_header_filename(filename,
9781     >                      hversion,hnfft,hunita,hispin,hne)
9782
9783        if ( (hnfft(1).ne.control_ngrid(1)) .or.
9784     >       (hnfft(2).ne.control_ngrid(2)) .or.
9785     >       (hnfft(3).ne.control_ngrid(3)) ) then
9786
9787        hnfft(1) = control_ngrid(1)
9788        hnfft(2) = control_ngrid(2)
9789        hnfft(3) = control_ngrid(3)
9790        call Parallel_taskid(taskid)
9791        value = btdb_parallel(.false.)
9792        if (taskid.eq.MASTER) then
9793
9794
9795          ind = index(filename,' ') - 1
9796          if (.not. btdb_cput(rtdb,'xpndr:old_wavefunction_filename',
9797     >                    1,filename(1:ind)))
9798     >     call errquit(
9799     >     'wvfnc_expander_input: btdb_cput failed', 0, RTDB_ERR)
9800
9801          if (.not. btdb_cput(rtdb,'xpndr:new_wavefunction_filename',
9802     >                    1,filename(1:ind)))
9803     >     call errquit(
9804     >     'wvfnc_expander_input: btdb_cput failed', 0, RTDB_ERR)
9805
9806          if (.not. btdb_put(rtdb,'xpndr:ngrid',mt_int,3,hnfft))
9807     >     call errquit(
9808     >     'wvfnc_expander_input: btdb_put failed', 0, RTDB_ERR)
9809
9810          write(*,*)
9811          write(*,*) "Grid is being converted:"
9812          write(*,*) "------------------------"
9813          write(*,*)
9814          write(*,*) "To turn off automatic grid conversion:"
9815          write(*,*)
9816          write(*,*) "set nwpw:psi_nogrid .false."
9817          write(*,*)
9818          value = wvfnc_expander(rtdb)
9819
9820        end if
9821        value = btdb_parallel(.true.)
9822
9823      end if
9824
9825      end if
9826
9827*     *****  read initial wavefunctions into psi1  ****
9828      call epsi_read(ispin0,ne_excited,dcpl_mb(psi1_excited(1)))
9829
9830
9831*     **** force inversion symmetry ****
9832      if (control_symmetry().eq.1)  then
9833         call Parallel_taskid(taskid)
9834         if (taskid.eq.MASTER) then
9835         write(*,*)
9836         write(*,*)
9837     >   "Projecting virtual wavefunctions to have inversion symmetry"
9838         write(*,*)
9839         end if
9840         call psi_symm_project(ispin0,ne_excited,npack1,
9841     >                         dcpl_mb(psi1_excited(1)))
9842      end if
9843
9844
9845c*     **** Ortho Check ****
9846c      call Grsm_gg_trace(npack1,(ne_excited(1)+ne_excited(2)),
9847c     >                        dcpl_mb(psi1_excited(1)),
9848c     >                        dcpl_mb(psi1_excited(1)),
9849c     >                        sum2)
9850c
9851c      sum1 = dble(ne_excited(1) + ne_excited(2))
9852c      if (dabs(sum2-sum1).gt.1.0d-10) then
9853c
9854c         call Parallel_taskid(taskid)
9855c         call Grsm_g_MakeOrtho(npack1,ne_excited(1),
9856c     >                         dcpl_mb(psi1_excited(1)))
9857c         if (ispin.gt.1) then
9858c           call Grsm_g_MakeOrtho(npack1,ne_excited(2),
9859c     >                           dcpl_mb(psi1_excited(1)
9860c     >                                  +ne_excited(1)*npack1))
9861c         end if
9862c         call Grsm_gg_trace(npack1,(ne_excited(1)+ne_excited(2)),
9863c     >                        dcpl_mb(psi1_excited(1)),
9864c     >                        dcpl_mb(psi1_excited(1)),
9865c     >                        sum2)
9866c         if (taskid.eq.MASTER)
9867c     >    write(*,*) "Warning - Gram-Schmidt being performed on epsi:",
9868c     >               dabs(sum2-sum1)
9869c
9870c      end if
9871
9872
9873      epsi_initialize = value
9874      return
9875      end
9876
9877
9878
9879*     ***************************
9880*     *				*
9881*     *		epsi_finalize	*
9882*     *				*
9883*     ***************************
9884
9885      logical function epsi_finalize(writepsi)
9886      implicit none
9887      logical writepsi
9888
9889#include "bafdecls.fh"
9890#include "errquit.fh"
9891#include "psi.fh"
9892
9893*     **** local variables ****
9894      logical value
9895
9896*     ***** write psi1 wavefunctions ****
9897      if (writepsi)
9898     >  call epsi_write(ispin,ne_excited,dcpl_mb(psi1_excited(1)))
9899
9900      value = BA_free_heap(hml_excited(2))
9901      value = value.and.BA_free_heap(eig_excited(2))
9902      value = value.and.BA_free_heap(psi2_excited(2))
9903      value = value.and.BA_free_heap(psi1_excited(2))
9904      if (.not. value)
9905     >  call errquit('epsi_finalize: error freeing heap',0, MA_ERR)
9906
9907      epsi_finalize = value
9908      return
9909      end
9910
9911
9912*     ***********************
9913*     *			    *
9914*     *	    psi_Mulliken    *
9915*     *			    *
9916*     ***********************
9917
9918      subroutine psi_Mulliken(rtdb)
9919      implicit none
9920      integer rtdb
9921
9922#include "bafdecls.fh"
9923#include "psi.fh"
9924#include "stdio.fh"
9925
9926
9927*     **** Lubin Water Analysis ****
9928      call pspw_Lubin_water_analysis(rtdb,ispin,ne,2*nfft3d,
9929     >                                 dbl_mb(rho1(1)))
9930
9931*     **** Atom Analysis ****
9932      call pspw_atom_analysis(rtdb,ispin,2*nfft3d,dbl_mb(rho1(1)))
9933
9934*     **** Mulliken Analysis ****
9935      call pspw_analysis(0,rtdb,ispin,ne,dcpl_mb(psi2(1)),
9936     >                                   dbl_mb(eig(1)))
9937
9938      call pspw_gen_APC(ispin,ne,dcpl_mb(dng1(1)),.false.)
9939      call pspw_print_APC(luout)
9940
9941      call pspw_gen_atom_Efield(rtdb,ispin,dbl_mb(rho1(1)),
9942     >                          dcpl_mb(dng1(1)))
9943
9944      call pspw_gen_atom_Efield_grad(rtdb,ispin,ne,
9945     >                          dcpl_mb(psi1(1)),
9946     >                          dcpl_mb(dng1(1)))
9947
9948      return
9949      end
9950
9951
9952*     ***********************
9953*     *                     *
9954*     *     psi_Born        *
9955*     *                     *
9956*     ***********************
9957
9958      subroutine psi_Born()
9959      implicit none
9960
9961#include "bafdecls.fh"
9962#include "psi.fh"
9963#include "stdio.fh"
9964
9965      call pspw_gen_APC(ispin,ne,dcpl_mb(dng1(1)),.false.)
9966      call pspw_print_APC(luout)
9967
9968      return
9969      end
9970
9971
9972
9973*     ***********************
9974*     *                     *
9975*     *    epsi_Mulliken    *
9976*     *                     *
9977*     ***********************
9978
9979      subroutine epsi_Mulliken(rtdb)
9980      implicit none
9981      integer rtdb
9982
9983#include "bafdecls.fh"
9984#include "psi.fh"
9985
9986      call pspw_analysis(1,rtdb,ispin,ne_excited,
9987     >                   dcpl_mb(psi1_excited(1)),
9988     >                   dbl_mb(eig_excited(1)))
9989      return
9990      end
9991
9992
9993
9994
9995*     ***********************
9996*     *                     *
9997*     *     psi_DOS         *
9998*     *                     *
9999*     ***********************
10000
10001      subroutine psi_DOS(rtdb)
10002      implicit none
10003      integer rtdb
10004
10005#include "btdb.fh"
10006#include "bafdecls.fh"
10007#include "psi.fh"
10008#include "errquit.fh"
10009
10010
10011*     **** local variables ****
10012      logical value
10013      integer npoints,ii
10014      integer weight(2),nemax
10015      real*8 emin,emax,alpha
10016      character*255 filename
10017
10018      nemax = ne(1)
10019      value = BA_push_get(mt_dbl,(nemax),'weight',weight(2),weight(1))
10020      if (.not. value)
10021     >  call errquit('psi_dos:out of stack memory',0, MA_ERR)
10022      call dcopy(nemax,1.0d0,0,dbl_mb(weight(1)),1)
10023
10024
10025      if (.not.btdb_get(rtdb,'dos:alpha',mt_dbl,1,alpha)) then
10026        alpha = 0.05d0/27.2116d0
10027      end if
10028
10029      if (.not.btdb_get(rtdb,'dos:npoints',mt_int,1,npoints)) then
10030        npoints = 500
10031      end if
10032
10033      if (.not.btdb_get(rtdb,'dos:emin',mt_dbl,1,emin)) then
10034         emin = 99999.0d0
10035         do ii=1,ne(1)+ne(2)
10036           if (dbl_mb(eig(1)+ii-1).lt.emin) emin = dbl_mb(eig(1)+ii-1)
10037         end do
10038         emin = emin - 0.1d0
10039      end if
10040
10041      if (.not.btdb_get(rtdb,'dos:emax',mt_dbl,1,emax)) then
10042         emax = -99999.0d0
10043         do ii=1,ne(1)+ne(2)
10044           if (dbl_mb(eig(1)+ii-1).gt.emax) emax = dbl_mb(eig(1)+ii-1)
10045         end do
10046         emax = emax + 0.1d0
10047      end if
10048
10049*     **** generate DENSITY OF STATES *****
10050      if (ispin.eq.1) then
10051        filename = "smear_dos_both"
10052        call densityofstates(filename,.false.,
10053     >                     dbl_mb(eig(1)),dbl_mb(weight(1)),ne(1),
10054     >                     1.0d0,alpha,npoints,emin,emax)
10055        filename = "smear_fdos_both"
10056        call densityofstates(filename,.false.,
10057     >                     dbl_mb(eig(1)),dbl_mb(weight(1)),ne(1),
10058     >                     1.0d0,alpha,npoints,emin,emax)
10059      end if
10060
10061      if (ispin.eq.2) then
10062        filename = "smear_dos_alpha"
10063        call densityofstates(filename,.false.,
10064     >                     dbl_mb(eig(1)),dbl_mb(weight(1)),ne(1),
10065     >                     1.0d0,alpha,npoints,emin,emax)
10066        filename = "smear_dos_beta"
10067        call densityofstates(filename,.false.,
10068     >               dbl_mb(eig(1)+ne(1)),dbl_mb(weight(1)),ne(2),
10069     >               -1.0d0,alpha,npoints,emin,emax)
10070        filename = "smear_fdos_alpha"
10071        call densityofstates(filename,.false.,
10072     >                     dbl_mb(eig(1)),dbl_mb(weight(1)),ne(1),
10073     >                     1.0d0,alpha,npoints,emin,emax)
10074        filename = "smear_fdos_beta"
10075        call densityofstates(filename,.false.,
10076     >               dbl_mb(eig(1)+ne(1)),dbl_mb(weight(1)),ne(2),
10077     >               -1.0d0,alpha,npoints,emin,emax)
10078      end if
10079
10080
10081      value = BA_pop_stack(weight(2))
10082      if (.not. value)
10083     >  call errquit('psi_dos: error freeing stack',0, MA_ERR)
10084
10085      return
10086      end
10087
10088*     ***********************
10089*     *                     *
10090*     *     epsi_DOS        *
10091*     *                     *
10092*     ***********************
10093
10094      subroutine epsi_DOS(rtdb)
10095      implicit none
10096      integer rtdb
10097
10098#include "btdb.fh"
10099#include "bafdecls.fh"
10100#include "psi.fh"
10101#include "errquit.fh"
10102
10103
10104*     **** local variables ****
10105      logical value
10106      integer weight(2),npoints,ii
10107      real*8 emin,emax,alpha
10108      character*255 filename
10109
10110      value = BA_push_get(mt_dbl,(ne_excited(1)+ne_excited(2)),
10111     >                    'weight',weight(2),weight(1))
10112      if (.not. value)
10113     >  call errquit('epsi_dos:out of stack memory',0, MA_ERR)
10114      call dcopy((ne_excited(1)+ne_excited(2)),
10115     >          1.0d0,0,dbl_mb(weight(1)),1)
10116
10117      if (.not.btdb_get(rtdb,'dos:alpha',mt_dbl,1,alpha)) then
10118        alpha = 0.05d0/27.2116d0
10119      end if
10120
10121      if (.not.btdb_get(rtdb,'dos:npoints',mt_int,1,npoints)) then
10122        npoints = 500
10123      end if
10124
10125      if (.not.btdb_get(rtdb,'dos:emin',mt_dbl,1,emin)) then
10126         emin = 99999.0d0
10127         do ii=1,ne_excited(1)+ne_excited(2)
10128           if (dbl_mb(eig_excited(1)+ii-1).lt.emin)
10129     >       emin = dbl_mb(eig_excited(1)+ii-1)
10130         end do
10131         emin = emin - 0.1d0
10132      end if
10133
10134      if (.not.btdb_get(rtdb,'dos:emax',mt_dbl,1,emax)) then
10135         emax = -99999.0d0
10136         do ii=1,ne_excited(1)+ne_excited(2)
10137           if (dbl_mb(eig_excited(1)+ii-1).gt.emax)
10138     >       emax = dbl_mb(eig_excited(1)+ii-1)
10139         end do
10140         emax = emax + 0.1d0
10141      end if
10142
10143*     **** generate DENSITY OF STATES *****
10144      if (ispin.eq.1) then
10145        filename = "smear_vdos_both"
10146        call densityofstates(filename,.false.,
10147     >                     dbl_mb(eig_excited(1)),dbl_mb(weight(1)),
10148     >                     ne_excited(1),
10149     >                     1.0d0,alpha,npoints,emin,emax)
10150        filename = "smear_dos_both"
10151        call densityofstates(filename,.true.,
10152     >                     dbl_mb(eig_excited(1)),dbl_mb(weight(1)),
10153     >                     ne_excited(1),
10154     >                     1.0d0,alpha,npoints,emin,emax)
10155      end if
10156
10157      if (ispin.eq.2) then
10158        filename = "smear_vdos_alpha"
10159        call densityofstates(filename,.false.,
10160     >                     dbl_mb(eig_excited(1)),dbl_mb(weight(1)),
10161     >                     ne_excited(1),
10162     >                     1.0d0,alpha,npoints,emin,emax)
10163        filename = "smear_vdos_beta"
10164        call densityofstates(filename,.false.,
10165     >           dbl_mb(eig_excited(1)+ne_excited(1)),dbl_mb(weight(1)),
10166     >           ne_excited(2),
10167     >           -1.0d0,alpha,npoints,emin,emax)
10168        filename = "smear_dos_alpha"
10169        call densityofstates(filename,.true.,
10170     >                     dbl_mb(eig_excited(1)),dbl_mb(weight(1)),
10171     >                     ne_excited(1),
10172     >                     1.0d0,alpha,npoints,emin,emax)
10173        filename = "smear_dos_beta"
10174        call densityofstates(filename,.true.,
10175     >           dbl_mb(eig_excited(1)+ne_excited(1)),dbl_mb(weight(1)),
10176     >           ne_excited(2),
10177     >           -1.0d0,alpha,npoints,emin,emax)
10178      end if
10179
10180      value = BA_pop_stack(weight(2))
10181      if (.not. value)
10182     >  call errquit('epsi_dos: error freeing stack',0, MA_ERR)
10183
10184      return
10185      end
10186
10187cccccccccccccccccccccccccccccccccccccccccccccccccccccc
10188      subroutine psi_polariz()
10189      implicit none
10190#include "psi.fh"
10191#include "bafdecls.fh"
10192#include "errquit.fh"
10193#include "util.fh"
10194      integer psirx(2),asize
10195      logical val
10196      asize=(ne(1)+ne(2))*nfft3d*2
10197      val=BA_push_get(mt_dbl,asize,"psir",psirx(2),psirx(1))
10198      if (.not.val) then
10199        call errquit("psi_polariz stack empty",0,0)
10200      end if
10201      call berry_phase_pol(psi2,psirx)
10202      val=BA_pop_stack(psirx(2))
10203      if (.not.val) then
10204        call errquit("psi_polariz failed pop stack",0,0)
10205      end if
10206      return
10207      end
10208ccccccccccccccccccccccccccccccccccccccccccccccccccc
10209
10210
10211      subroutine epsi_generate_kb_vnm(vnm)
10212      implicit none
10213      real*8 vnm(*)
10214
10215#include "bafdecls.fh"
10216#include "errquit.fh"
10217#include "psi.fh"
10218
10219*     **** local variables ****
10220      logical value
10221      integer ms,i,j,vshift,ishift,neall(2)
10222      integer Gx(2),Gy(2),Gz(2),tmp(2),psii_ptr,psij_ptr
10223      real*8  vv(3)
10224
10225*     **** external functions ****
10226      integer  Dneall_mne_size,G_indx
10227      external Dneall_mne_size,G_indx
10228
10229      neall(1) = ne(1)+ne_excited(1)
10230      neall(2) = ne(2)+ne_excited(2)
10231
10232      value = BA_push_get(mt_dbl, nfft3d,'Gx',Gx(2),Gx(1))
10233      value = value.and.
10234     >        BA_push_get(mt_dbl, nfft3d,'Gy',Gy(2),Gy(1))
10235      value = value.and.
10236     >        BA_push_get(mt_dbl, nfft3d,'Gz',Gz(2),Gz(1))
10237      value = value.and.
10238     >        BA_push_get(mt_dbl, 2*npack1,'tmp',tmp(2),tmp(1))
10239         if (.not. value)
10240     >      call errquit('epsi_generate_kb_vnm:pushing stack',1,MA_ERR)
10241
10242
10243*     **** define Gx,Gy and Gz in packed space ****
10244      call D3dB_t_Copy(1,dbl_mb(G_indx(1)),dbl_mb(Gx(1)))
10245      call D3dB_t_Copy(1,dbl_mb(G_indx(2)),dbl_mb(Gy(1)))
10246      call D3dB_t_Copy(1,dbl_mb(G_indx(3)),dbl_mb(Gz(1)))
10247      call Pack_t_pack(1,dbl_mb(Gx(1)))
10248      call Pack_t_pack(1,dbl_mb(Gy(1)))
10249      call Pack_t_pack(1,dbl_mb(Gz(1)))
10250
10251
10252      vshift = Dneall_mne_size(0,neall)
10253      !call dcopy(3*vshift,0.0d0,0,vnm,1)
10254      call Parallel_shared_vector_zero(.true.,3*vshift,vnm)
10255
10256      do ms=1,ispin
10257         ishift = (ms-1)*ne(1)
10258
10259         do j=1,neall(ms)
10260            if (j.le.ne(ms)) then
10261               psij_ptr =(j-1+ishift)*npack1 + psi1(1)
10262            else
10263               psij_ptr =(j-ne(ms)-1+ishift)*npack1 + psi1_excited(1)
10264            end if
10265
10266            do i=j,neall(ms)
10267               if (i.le.ne(ms)) then
10268                  psii_ptr =(i-1+ishift)*npack1 + psi1(1)
10269               else
10270                  psii_ptr =(i-ne(ms)-1+ishift)*npack1 + psi1_excited(1)
10271               end if
10272
10273               call Pack_tc_Mul(1,
10274     >                          dbl_mb(Gx(1)),
10275     >                          dcpl_mb(psij_ptr),
10276     >                          dbl_mb(tmp(1)))
10277               call Pack_cc_dot(1,
10278     >                          dcpl_mb(psii_ptr),
10279     >                          dbl_mb(tmp(1)),
10280     >                          vv(1))
10281               call Pack_tc_Mul(1,
10282     >                          dbl_mb(Gy(1)),
10283     >                          dcpl_mb(psij_ptr),
10284     >                          dbl_mb(tmp(1)))
10285               call Pack_cc_dot(1,
10286     >                          dcpl_mb(psii_ptr),
10287     >                          dbl_mb(tmp(1)),
10288     >                          vv(2))
10289               call Pack_tc_Mul(1,
10290     >                          dbl_mb(Gz(1)),
10291     >                          dcpl_mb(psij_ptr),
10292     >                          dbl_mb(tmp(1)))
10293               call Pack_cc_dot(1,
10294     >                          dcpl_mb(psii_ptr),
10295     >                          dbl_mb(tmp(1)),
10296     >                          vv(3))
10297
10298
10299               call Dneall_mne_set_value(vv(1),0,neall,ms,i,j,vnm)
10300               call Dneall_mne_set_value(vv(2),0,neall,ms,i,j,
10301     >                                   vnm(1+vshift))
10302               call Dneall_mne_set_value(vv(3),0,neall,ms,i,j,
10303     >                                   vnm(1+vshift+vshift))
10304               if (i.ne.j) then
10305                  call Dneall_mne_set_value(vv(1),0,neall,ms,j,i,vnm)
10306                  call Dneall_mne_set_value(vv(2),0,neall,ms,j,i,
10307     >                                      vnm(1+vshift))
10308                  call Dneall_mne_set_value(vv(3),0,neall,ms,j,i,
10309     >                                      vnm(1+vshift+vshift))
10310               end if
10311            end do
10312         end do
10313      end do
10314
10315      value =           BA_pop_stack(tmp(2))
10316      value = value.and.BA_pop_stack(Gz(2))
10317      value = value.and.BA_pop_stack(Gy(2))
10318      value = value.and.BA_pop_stack(Gx(2))
10319      if (.not. value)
10320     >   call errquit('epsi_generate_kb_vnm:popping stack',1,MA_ERR)
10321
10322
10323      return
10324      end
10325
10326
10327
10328*     *********************************
10329*     *                               *
10330*     *     psi_1pressure_stress      *
10331*     *                               *
10332*     *********************************
10333
10334      subroutine psi_1pressure_stress(pressure,p1,p2,stress)
10335      implicit none
10336      real*8 pressure,p1,p2,stress(3,3)
10337
10338#include "bafdecls.fh"
10339#include "errquit.fh"
10340#include "psi.fh"
10341
10342*     ***** rhoall common block ****
10343      integer rho1_all(2)
10344      integer rho2_all(2)
10345      common / rhoall_block / rho1_all,rho2_all
10346
10347*     ***** external functions *****
10348      integer  electron_xcp_ptr
10349      external electron_xcp_ptr
10350      real*8   psi_1vnl,rho_1exc,rho_1pxc
10351      external psi_1vnl,rho_1exc,rho_1pxc
10352
10353      call cgsd_pressure_stress(ispin,neq,
10354     >                          dcpl_mb(psi1(1)),
10355     >                          dbl_mb(rho1_all(1)),
10356     >                          dcpl_mb(dng1(1)),
10357     >                          dbl_mb(electron_xcp_ptr()),
10358     >                          psi_1vnl(),rho_1exc(),rho_1pxc(),
10359     >                          pressure,p1,p2,stress)
10360
10361
10362      return
10363      end
10364
10365
10366
10367*     ***********************
10368*     *                     *
10369*     *    psi_MP2_energy   *
10370*     *                     *
10371*     ***********************
10372      subroutine psi_MP2_energy(rtdb)
10373      implicit none
10374      integer rtdb
10375
10376#include "stdio.fh"
10377#include "bafdecls.fh"
10378#include "errquit.fh"
10379#include "psi.fh"
10380
10381*     **** local variables ****
10382      integer taskid,MASTER
10383      parameter (MASTER=0)
10384
10385      logical oprint,value
10386      integer ms,a,b,r,s,n2ft3d,icount
10387      real*8 ea,eb,er,es
10388      real*8 e2,d2,tmp2
10389      integer vpsi_r(2),v1h(2),v2h(2)
10390
10391*     **** allocate memory from heap ****
10392      n2ft3d = 2*nfft3d
10393      value = BA_alloc_get(mt_dcpl,nfft3d,'v1h',v1h(2),v1h(1))
10394      value = value.and.
10395     >        BA_alloc_get(mt_dcpl,nfft3d,'v2h',v2h(2),v2h(1))
10396      value = value.and.
10397     >        BA_alloc_get(mt_dbl,n2ft3d*(ne_excited(1)+ne_excited(2)),
10398     >                     'vpsi_r',vpsi_r(2),vpsi_r(1))
10399      if (.not. value)
10400     >  call errquit('psi_MP2_energy: error allocating heap memory',0,
10401     &       MA_ERR)
10402
10403
10404      call Parallel_taskid(taskid)
10405      oprint = (taskid.eq.MASTER)
10406
10407      icount = 0
10408      do ms=1,ispin
10409         do a=1,ne(ms)-1
10410            ea = dbl_mb(eig(1)+a-1+(ms-1)*ne(1))
10411            do r=1,ne_excited(ms)-1
10412               er = dbl_mb(eig_excited(1)+r-1+(ms-1)*ne(1))
10413
10414               do b=a+1,ne(ms)
10415                  eb = dbl_mb(eig(1)+b-1+(ms-1)*ne(1))
10416                  do s=r+1,ne_excited(ms)
10417                     es = dbl_mb(eig_excited(1)+s-1+(ms-1)*ne(1))
10418                     d2 = ea+eb-er-es
10419                     tmp2 = 1.0d0
10420                     e2 = e2 + tmp2*tmp2/d2
10421                     icount = icount + 1
10422                  if (oprint)
10423     >            write(luout,*) "a,b,r,s, Esub=",a,b,r,s,(tmp2*tmp2/d2)
10424                  end do
10425               end do
10426            end do
10427         end do
10428      end do
10429      if (ispin.eq.1) e2=e2+e2
10430
10431*     **** deallocate memory from heap ****
10432      value =     BA_free_heap(v1h(2))
10433     >       .and.BA_free_heap(v2h(2))
10434     >       .and.BA_free_heap(vpsi_r(2))
10435      if (.not. value)
10436     >  call errquit('psi_MP2_energy: error freeing heap memory',1,
10437     &       MA_ERR)
10438
10439
10440
10441      if (oprint) then
10442         write(luout,*) "EMP2 = ", e2, " icount=",icount
10443      end if
10444
10445      return
10446      end
10447
10448
10449
10450*     ***********************
10451*     *                     *
10452*     *  psi_2q_integrals   *
10453*     *                     *
10454*     ***********************
10455      subroutine psi_2q_integrals(rtdb)
10456      implicit none
10457      integer rtdb
10458
10459#include "stdio.fh"
10460#include "bafdecls.fh"
10461#include "errquit.fh"
10462#include "psi.fh"
10463
10464*     **** local variables ****
10465      integer taskid,MASTER,taskid_j,pj
10466      parameter (MASTER=0)
10467
10468      logical oprint,value
10469      integer icount,iicount,version
10470      integer i,j,k,l,ij,kl
10471      integer nall,nnall,n4all,nx,ny,nz
10472      integer psiij
10473      integer ipackm(2),jpackm(2)
10474      integer h1_integrals(2),h2_integrals(2)
10475      integer vij(2),dnij(2),orbi(2),orbj(2)
10476      real*8 e1,e2,scal1,scal2,dv
10477
10478*     **** external functions ****
10479      integer  control_version
10480      external control_version
10481      real*8   lattice_omega
10482      external lattice_omega
10483
10484      call Parallel_taskid(taskid)
10485      call Parallel2d_taskid_j(taskid_j)
10486      oprint = (taskid.eq.MASTER)
10487      version = control_version()
10488      if (version.eq.4) then
10489
10490         if (oprint) then
10491            write(luout,*)
10492     >      "== Generating One-Electron and Two-Electron Integrals =="
10493            write(luout,*)
10494         end if
10495         call D3dB_nx(1,nx)
10496         call D3dB_ny(1,ny)
10497         call D3dB_nz(1,nz)
10498         scal1 = 1.0d0/dble(nx*ny*nz)
10499         scal2 = 1.0d0/lattice_omega()
10500         dv    = scal1*lattice_omega()
10501
10502         nall = ne(1) + ne_excited(1)
10503         nnall = (nall*(nall+1))/2
10504         n4all = (nnall*(nnall+1))/2
10505
10506         value = BA_alloc_get(mt_dcpl,nfft3d,'vij',vij(2),vij(1))
10507         value = value.and.
10508     >           BA_alloc_get(mt_dcpl,nfft3d,'dnij',dnij(2),dnij(1))
10509         value = value.and.
10510     >           BA_alloc_get(mt_dcpl,nfft3d,'orbi',orbi(2),orbi(1))
10511         value = value.and.
10512     >           BA_alloc_get(mt_dcpl,nfft3d,'orbj',orbj(2),orbj(1))
10513         value = value.and.
10514     >           BA_alloc_get(mt_dbl,nnall,'h1_integrals',
10515     >                        h1_integrals(2),h1_integrals(1))
10516         value = value.and.
10517     >           BA_alloc_get(mt_dbl,n4all,'h2_integrals',
10518     >                        h2_integrals(2),h2_integrals(1))
10519         value = value.and.
10520     >           BA_alloc_get(mt_int,nnall,'ipackm',
10521     >                        ipackm(2),ipackm(1))
10522         value = value.and.
10523     >           BA_alloc_get(mt_int,nnall,'jpackm',
10524     >                        jpackm(2),jpackm(1))
10525         if (.not.value)
10526     >   call errquit('psi_2q_integrals:error allocating heap',0,MA_ERR)
10527
10528         call Parallel_shared_vector_zero(.true.,nnall,
10529     >                                    dbl_mb(h1_integrals(1)))
10530         call Parallel_shared_vector_zero(.true.,n4all,
10531     >                                    dbl_mb(h2_integrals(1)))
10532
10533         !**** generate 1e integrals - H1 only contains kinetic + ion-electron potentials ****
10534         call electron_gen_vall0()
10535         icount = 0
10536         do i=1,nall
10537            call psi_fetch_orbi_replicated(i,dcpl_mb(orbi(1)))
10538            call electron_get_H0psi_k_orb(dcpl_mb(orbi(1)),
10539     >                                    dcpl_mb(vij(1)))
10540
10541            do j=i,nall
10542               call psi_fetch_orbj_indx(j,psiij,pj)
10543               if (pj.eq.taskid_j) then
10544                  call Pack_cc_idot(1,dcpl_mb(psiij),dcpl_mb(vij(1)),e1)
10545                  dbl_mb(h1_integrals(1)+icount) = e1
10546               end if
10547               int_mb(ipackm(1) + icount) = i
10548               int_mb(jpackm(1) + icount) = j
10549               icount = icount + 1
10550            end do
10551            call ga_sync()
10552         end do
10553         call Parallel_Vector_SumAll(nnall,dbl_mb(h1_integrals(1)))
10554
10555         if (oprint) then
10556            write(luout,*) "begin_one_electron_integrals"
10557            icount = 0
10558            do i=1,nall
10559               do j=i,nall
10560                  write(luout,'(I5,I5,F21.10)')
10561     >            j,i,dbl_mb(h1_integrals(1)+icount)
10562                  icount = icount + 1
10563               end do
10564            end do
10565            write(luout,*) "end_one_electron_integrals"
10566         end if
10567
10568         !*** generate 2e integrals ***
10569         iicount = 0
10570         do ij=1,nnall
10571            i = int_mb(ipackm(1)+ij-1)
10572            j = int_mb(jpackm(1)+ij-1)
10573            call psi_fetch_orbi_replicated(i,dcpl_mb(orbi(1)))
10574            call Pack_c_unpack(1,   dcpl_mb(orbi(1)))
10575            call D3dB_cr_pfft3b(1,1,dcpl_mb(orbi(1)))
10576            call psi_fetch_orbi_replicated(j,dcpl_mb(orbj(1)))
10577            call Pack_c_unpack(1,   dcpl_mb(orbj(1)))
10578            call D3dB_cr_pfft3b(1,1,dcpl_mb(orbj(1)))
10579
10580*           **** generate dnij for Vij  ****
10581            call D3dB_rr_Mul(1,dcpl_mb(orbi(1)),
10582     >                         dcpl_mb(orbj(1)),
10583     >                         dcpl_mb(dnij(1)))
10584            call D3dB_r_SMul1(1,scal2,dcpl_mb(dnij(1)))
10585            call coulomb2_v(dcpl_mb(dnij(1)),dcpl_mb(vij(1)))
10586
10587
10588            do kl=ij,nnall
10589               k = int_mb(ipackm(1)+kl-1)
10590               l = int_mb(jpackm(1)+kl-1)
10591               call psi_fetch_orbi_replicated(k,dcpl_mb(orbi(1)))
10592               call Pack_c_unpack(1,   dcpl_mb(orbi(1)))
10593               call D3dB_cr_pfft3b(1,1,dcpl_mb(orbi(1)))
10594               call psi_fetch_orbi_replicated(l,dcpl_mb(orbj(1)))
10595               call Pack_c_unpack(1,   dcpl_mb(orbj(1)))
10596               call D3dB_cr_pfft3b(1,1,dcpl_mb(orbj(1)))
10597
10598               call D3dB_rr_Mul(1,dcpl_mb(orbi(1)),
10599     >                            dcpl_mb(orbj(1)),
10600     >                            dcpl_mb(dnij(1)))
10601               call D3dB_r_SMul1(1,scal2,dcpl_mb(dnij(1)))
10602
10603               call D3dB_rr_dot(1,dcpl_mb(dnij(1)),dcpl_mb(vij(1)),e2)
10604               !e2 = 0.5d0*e2*dv
10605               e2 = e2*dv
10606               dbl_mb(h2_integrals(1)+iicount) = e2
10607
10608               !write(*,*) "i,j,k,l,e2=",i,j,k,l,e2
10609               iicount = iicount + 1
10610            end do
10611         end do
10612
10613
10614         if (oprint) then
10615            write(luout,*)
10616            write(luout,*) "begin_two_electron_integrals"
10617            iicount = 0
10618            do ij=1,nnall
10619               i = int_mb(ipackm(1)+ij-1)
10620               j = int_mb(jpackm(1)+ij-1)
10621               do kl=ij,nnall
10622                  k = int_mb(ipackm(1)+kl-1)
10623                  l = int_mb(jpackm(1)+kl-1)
10624                  write(luout,'(I5,I5,I5,I5,F20.10)')
10625     >            j,i,l,k,dbl_mb(h2_integrals(1)+iicount)
10626                  iicount = iicount + 1
10627               end do
10628            end do
10629            write(luout,*) "end_two_electron_integrals"
10630            write(luout,*)
10631         end if
10632
10633
10634*        **** deallocate memory from heap ****
10635         value = BA_free_heap(vij(2))
10636         value = value.and.BA_free_heap(dnij(2))
10637         value = value.and.BA_free_heap(orbi(2))
10638         value = value.and.BA_free_heap(orbj(2))
10639         value = value.and.BA_free_heap(h1_integrals(2))
10640         value = value.and.BA_free_heap(h2_integrals(2))
10641         value = value.and.BA_free_heap(ipackm(2))
10642         value = value.and.BA_free_heap(jpackm(2))
10643         if (.not. value)
10644     >   call errquit('psi_2q_integrals:error freeing heap',1,MA_ERR)
10645      end if
10646
10647      return
10648      end
10649
10650      subroutine psi_fetch_orbi_replicated(i,orbi)
10651      implicit none
10652      integer i
10653      complex*16 orbi(*)
10654
10655#include "bafdecls.fh"
10656#include "errquit.fh"
10657#include "psi.fh"
10658
10659*     **** local variables ****
10660      integer pi,q,taskid_j,psii_ptr
10661
10662      call Parallel2d_taskid_j(taskid_j)
10663
10664      if (i.le.ne(1)) then
10665         call Dneall_ntoqp(i,q,pi)
10666         call Pack_c_Zero(1,orbi)
10667         if (pi.eq.taskid_j) then
10668            psii_ptr = psi1(1) + (q-1)*npack1
10669            call Pack_c_Copy(1,dcpl_mb(psii_ptr),orbi)
10670         end if
10671         call D1dB_Vector_SumAll(2*npack1,orbi)
10672      else
10673         psii_ptr = psi1_excited(1) + (i-1-ne(1))*npack1
10674         call Pack_c_Copy(1,dcpl_mb(psii_ptr),orbi)
10675      end if
10676
10677      return
10678      end
10679
10680      subroutine psi_fetch_orbj_indx(j,orbj_indx,pj)
10681      implicit none
10682      integer j
10683      integer orbj_indx
10684      integer pj
10685
10686#include "bafdecls.fh"
10687#include "errquit.fh"
10688#include "psi.fh"
10689
10690*     **** local variables ****
10691      integer q,np_j
10692
10693      call Parallel2d_np_j(np_j)
10694      if (j.le.ne(1)) then
10695         call Dneall_ntoqp(j,q,pj)
10696         orbj_indx = psi1(1) + (q-1)*npack1
10697      else
10698         orbj_indx = psi1_excited(1) + (j-1-ne(1))*npack1
10699         pj = mod(pj+1,np_j)
10700      end if
10701
10702      return
10703      end
10704
10705
10706