-
Notifications
You must be signed in to change notification settings - Fork 1
/
qs_initial_guess.F
1469 lines (1306 loc) · 67.4 KB
/
qs_initial_guess.F
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!--------------------------------------------------------------------------------------------------!
! CP2K: A general program to perform molecular dynamics simulations !
! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
! !
! SPDX-License-Identifier: GPL-2.0-or-later !
!--------------------------------------------------------------------------------------------------!
! **************************************************************************************************
!> \brief Routines to somehow generate an initial guess
!> \par History
!> 2006.03 Moved here from qs_scf.F [Joost VandeVondele]
! **************************************************************************************************
MODULE qs_initial_guess
USE atom_kind_orbitals, ONLY: calculate_atomic_orbitals
USE atomic_kind_types, ONLY: atomic_kind_type,&
get_atomic_kind,&
get_atomic_kind_set
USE basis_set_types, ONLY: get_gto_basis_set,&
gto_basis_set_type
USE cp_control_types, ONLY: dft_control_type
USE cp_dbcsr_api, ONLY: &
dbcsr_checksum, dbcsr_copy, dbcsr_dot, dbcsr_filter, dbcsr_get_diag, dbcsr_get_num_blocks, &
dbcsr_get_occupation, dbcsr_iterator_blocks_left, dbcsr_iterator_next_block, &
dbcsr_iterator_start, dbcsr_iterator_stop, dbcsr_iterator_type, dbcsr_multiply, &
dbcsr_nfullrows_total, dbcsr_p_type, dbcsr_release, dbcsr_scale, dbcsr_set, &
dbcsr_set_diag, dbcsr_type, dbcsr_verify_matrix
USE cp_dbcsr_operations, ONLY: copy_dbcsr_to_fm,&
copy_fm_to_dbcsr,&
cp_dbcsr_sm_fm_multiply,&
cp_fm_to_dbcsr_row_template
USE cp_fm_cholesky, ONLY: cp_fm_cholesky_decompose
USE cp_fm_struct, ONLY: cp_fm_struct_create,&
cp_fm_struct_get,&
cp_fm_struct_release,&
cp_fm_struct_type
USE cp_fm_types, ONLY: &
cp_fm_create, cp_fm_get_info, cp_fm_get_submatrix, cp_fm_init_random, cp_fm_release, &
cp_fm_set_all, cp_fm_set_submatrix, cp_fm_to_fm, cp_fm_type
USE cp_log_handling, ONLY: cp_get_default_logger,&
cp_logger_get_default_io_unit,&
cp_logger_type,&
cp_to_string
USE cp_output_handling, ONLY: cp_print_key_finished_output,&
cp_print_key_unit_nr
USE external_potential_types, ONLY: all_potential_type,&
gth_potential_type,&
sgp_potential_type
USE hfx_types, ONLY: hfx_type
USE input_constants, ONLY: &
atomic_guess, core_guess, eht_guess, history_guess, mopac_guess, no_guess, random_guess, &
restart_guess, sparse_guess
USE input_cp2k_hfx, ONLY: ri_mo
USE input_section_types, ONLY: section_vals_get_subs_vals,&
section_vals_type,&
section_vals_val_get
USE kinds, ONLY: default_path_length,&
dp
USE kpoint_io, ONLY: read_kpoints_restart
USE kpoint_types, ONLY: kpoint_type
USE message_passing, ONLY: mp_para_env_type
USE particle_methods, ONLY: get_particle_set
USE particle_types, ONLY: particle_type
USE qs_atomic_block, ONLY: calculate_atomic_block_dm
USE qs_density_matrices, ONLY: calculate_density_matrix
USE qs_dftb_utils, ONLY: get_dftb_atom_param
USE qs_eht_guess, ONLY: calculate_eht_guess
USE qs_environment_types, ONLY: get_qs_env,&
qs_environment_type
USE qs_kind_types, ONLY: get_qs_kind,&
get_qs_kind_set,&
qs_kind_type
USE qs_mo_io, ONLY: read_mo_set_from_restart,&
wfn_restart_file_name
USE qs_mo_methods, ONLY: make_basis_lowdin,&
make_basis_simple,&
make_basis_sm
USE qs_mo_occupation, ONLY: set_mo_occupation
USE qs_mo_types, ONLY: get_mo_set,&
mo_set_restrict,&
mo_set_type,&
reassign_allocated_mos
USE qs_mom_methods, ONLY: do_mom_guess
USE qs_rho_methods, ONLY: qs_rho_update_rho
USE qs_rho_types, ONLY: qs_rho_get,&
qs_rho_type
USE qs_scf_methods, ONLY: eigensolver,&
eigensolver_simple
USE qs_scf_types, ONLY: block_davidson_diag_method_nr,&
block_krylov_diag_method_nr,&
general_diag_method_nr,&
ot_diag_method_nr,&
qs_scf_env_type
USE qs_wf_history_methods, ONLY: wfi_update
USE scf_control_types, ONLY: scf_control_type
USE util, ONLY: sort
USE xtb_types, ONLY: get_xtb_atom_param,&
xtb_atom_type
#include "./base/base_uses.f90"
IMPLICIT NONE
PRIVATE
CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_initial_guess'
PUBLIC :: calculate_first_density_matrix, calculate_mopac_dm
PUBLIC :: calculate_atomic_fock_matrix
TYPE atom_matrix_type
REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: mat => NULL()
END TYPE atom_matrix_type
CONTAINS
! **************************************************************************************************
!> \brief can use a variety of methods to come up with an initial
!> density matrix and optionally an initial wavefunction
!> \param scf_env SCF environment information
!> \param qs_env QS environment
!> \par History
!> 03.2006 moved here from qs_scf [Joost VandeVondele]
!> 06.2007 allow to skip the initial guess [jgh]
!> 08.2014 kpoints [JGH]
!> 10.2019 tot_corr_zeff, switch_surf_dip [SGh]
!> \note
!> badly needs to be split in subroutines each doing one of the possible
!> schemes
! **************************************************************************************************
SUBROUTINE calculate_first_density_matrix(scf_env, qs_env)
TYPE(qs_scf_env_type), POINTER :: scf_env
TYPE(qs_environment_type), POINTER :: qs_env
CHARACTER(LEN=*), PARAMETER :: routineN = 'calculate_first_density_matrix'
CHARACTER(LEN=default_path_length) :: file_name, filename
INTEGER :: atom_a, blk, density_guess, handle, homo, i, iatom, ic, icol, id_nr, ikind, irow, &
iseed(4), ispin, istart_col, istart_row, j, last_read, n, n_cols, n_rows, nao, natom, &
natoms, natoms_tmp, nblocks, nelectron, nmo, nmo_tmp, not_read, nsgf, nspin, nvec, ounit, &
safe_density_guess, size_atomic_kind_set, z
INTEGER, ALLOCATABLE, DIMENSION(:) :: first_sgf, kind_of, last_sgf
INTEGER, DIMENSION(2) :: nelectron_spin
INTEGER, DIMENSION(:), POINTER :: atom_list, elec_conf, nelec_kind, &
sort_kind
LOGICAL :: did_guess, do_hfx_ri_mo, do_kpoints, do_std_diag, exist, has_unit_metric, &
natom_mismatch, need_mos, need_wm, ofgpw, owns_ortho, print_history_log, print_log
REAL(dp), ALLOCATABLE, DIMENSION(:, :) :: buff, buff2
REAL(dp), DIMENSION(:, :), POINTER :: pdata
REAL(KIND=dp) :: checksum, eps, length, maxocc, occ, &
rscale, tot_corr_zeff, trps1, zeff
REAL(KIND=dp), DIMENSION(0:3) :: edftb
TYPE(atom_matrix_type), DIMENSION(:), POINTER :: pmat
TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
TYPE(atomic_kind_type), POINTER :: atomic_kind
TYPE(cp_fm_struct_type), POINTER :: ao_ao_struct, ao_mo_struct
TYPE(cp_fm_type) :: sv
TYPE(cp_fm_type), DIMENSION(:), POINTER :: work1
TYPE(cp_fm_type), POINTER :: mo_coeff, moa, mob, ortho, work2
TYPE(cp_logger_type), POINTER :: logger
TYPE(dbcsr_iterator_type) :: iter
TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: h_core_sparse, matrix_ks, p_rmpv, &
s_sparse
TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_h_kp, matrix_ks_kp, matrix_s_kp, &
rho_ao_kp
TYPE(dbcsr_type) :: mo_dbcsr, mo_tmp_dbcsr
TYPE(dft_control_type), POINTER :: dft_control
TYPE(gto_basis_set_type), POINTER :: orb_basis_set
TYPE(hfx_type), DIMENSION(:, :), POINTER :: x_data
TYPE(kpoint_type), POINTER :: kpoints
TYPE(mo_set_type), DIMENSION(:), POINTER :: mo_array, mos_last_converged
TYPE(mp_para_env_type), POINTER :: para_env
TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
TYPE(qs_kind_type), POINTER :: qs_kind
TYPE(qs_rho_type), POINTER :: rho
TYPE(scf_control_type), POINTER :: scf_control
TYPE(section_vals_type), POINTER :: dft_section, input, subsys_section
logger => cp_get_default_logger()
NULLIFY (atomic_kind, qs_kind, mo_coeff, orb_basis_set, atomic_kind_set, &
qs_kind_set, particle_set, ortho, work2, work1, mo_array, s_sparse, &
scf_control, dft_control, p_rmpv, para_env, h_core_sparse, matrix_ks, rho, &
mos_last_converged)
NULLIFY (dft_section, input, subsys_section)
NULLIFY (matrix_s_kp, matrix_h_kp, matrix_ks_kp, rho_ao_kp)
NULLIFY (moa, mob)
NULLIFY (atom_list, elec_conf, kpoints)
edftb = 0.0_dp
tot_corr_zeff = 0.0_dp
CALL timeset(routineN, handle)
CALL get_qs_env(qs_env, &
atomic_kind_set=atomic_kind_set, &
qs_kind_set=qs_kind_set, &
particle_set=particle_set, &
mos=mo_array, &
matrix_s_kp=matrix_s_kp, &
matrix_h_kp=matrix_h_kp, &
matrix_ks_kp=matrix_ks_kp, &
input=input, &
scf_control=scf_control, &
dft_control=dft_control, &
has_unit_metric=has_unit_metric, &
do_kpoints=do_kpoints, &
kpoints=kpoints, &
rho=rho, &
nelectron_spin=nelectron_spin, &
para_env=para_env, &
x_data=x_data)
CALL qs_rho_get(rho, rho_ao_kp=rho_ao_kp)
IF (dft_control%switch_surf_dip) THEN
CALL get_qs_env(qs_env, mos_last_converged=mos_last_converged)
END IF
! just initialize the first image, the other density are set to zero
DO ispin = 1, dft_control%nspins
DO ic = 1, SIZE(rho_ao_kp, 2)
CALL dbcsr_set(rho_ao_kp(ispin, ic)%matrix, 0.0_dp)
END DO
END DO
s_sparse => matrix_s_kp(:, 1)
h_core_sparse => matrix_h_kp(:, 1)
matrix_ks => matrix_ks_kp(:, 1)
p_rmpv => rho_ao_kp(:, 1)
work1 => scf_env%scf_work1
work2 => scf_env%scf_work2
ortho => scf_env%ortho
dft_section => section_vals_get_subs_vals(input, "DFT")
nspin = dft_control%nspins
ofgpw = dft_control%qs_control%ofgpw
density_guess = scf_control%density_guess
do_std_diag = .FALSE.
do_hfx_ri_mo = .FALSE.
IF (ASSOCIATED(x_data)) THEN
IF (x_data(1, 1)%do_hfx_ri) THEN
IF (x_data(1, 1)%ri_data%flavor == ri_mo) do_hfx_ri_mo = .TRUE.
END IF
END IF
IF (ASSOCIATED(scf_env%krylov_space)) do_std_diag = (scf_env%krylov_space%eps_std_diag > 0.0_dp)
need_mos = scf_control%use_ot .OR. scf_env%method == ot_diag_method_nr .OR. &
(scf_env%method == block_krylov_diag_method_nr .AND. .NOT. do_std_diag) &
.OR. dft_control%do_admm .OR. scf_env%method == block_davidson_diag_method_nr &
.OR. do_hfx_ri_mo
safe_density_guess = atomic_guess
IF (dft_control%qs_control%semi_empirical .OR. dft_control%qs_control%dftb) THEN
IF (density_guess == atomic_guess) density_guess = mopac_guess
safe_density_guess = mopac_guess
END IF
IF (dft_control%qs_control%xtb) THEN
IF (do_kpoints) THEN
IF (density_guess == atomic_guess) density_guess = mopac_guess
safe_density_guess = mopac_guess
ELSE
IF (density_guess == atomic_guess) density_guess = core_guess
safe_density_guess = core_guess
END IF
END IF
IF (scf_control%use_ot .AND. &
(.NOT. ((density_guess == random_guess) .OR. &
(density_guess == atomic_guess) .OR. &
(density_guess == core_guess) .OR. &
(density_guess == mopac_guess) .OR. &
(density_guess == eht_guess) .OR. &
(density_guess == sparse_guess) .OR. &
(((density_guess == restart_guess) .OR. &
(density_guess == history_guess)) .AND. &
(scf_control%level_shift == 0.0_dp))))) THEN
CALL cp_abort(__LOCATION__, &
"OT needs GUESS ATOMIC / CORE / RANDOM / SPARSE / RESTART / HISTORY RESTART: other options NYI")
END IF
! if a restart was requested, check that the file exists,
! if not we fall back to an atomic guess. No kidding, the file name should remain
! in sync with read_mo_set_from_restart
id_nr = 0
IF (density_guess == restart_guess) THEN
! only check existence on I/O node, otherwise if file exists there but
! not on compute nodes, everything goes crazy even though only I/O
! node actually reads the file
IF (do_kpoints) THEN
IF (para_env%is_source()) THEN
CALL wfn_restart_file_name(file_name, exist, dft_section, logger, kp=.TRUE.)
END IF
ELSE
IF (para_env%is_source()) THEN
CALL wfn_restart_file_name(file_name, exist, dft_section, logger)
END IF
END IF
CALL para_env%bcast(exist)
CALL para_env%bcast(file_name)
IF (.NOT. exist) THEN
CALL cp_warn(__LOCATION__, &
"User requested to restart the wavefunction from the file named: "// &
TRIM(file_name)//". This file does not exist. Please check the existence of"// &
" the file or change properly the value of the keyword WFN_RESTART_FILE_NAME."// &
" Calculation continues using ATOMIC GUESS. ")
density_guess = safe_density_guess
END IF
ELSE IF (density_guess == history_guess) THEN
IF (do_kpoints) THEN
CPABORT("calculate_first_density_matrix: history_guess not implemented for k-points")
END IF
IF (para_env%is_source()) THEN
CALL wfn_restart_file_name(file_name, exist, dft_section, logger)
END IF
CALL para_env%bcast(exist)
CALL para_env%bcast(file_name)
nvec = qs_env%wf_history%memory_depth
not_read = nvec + 1
! At this level we read the saved backup RESTART files..
DO i = 1, nvec
j = i - 1
filename = TRIM(file_name)
IF (j /= 0) THEN
filename = TRIM(file_name)//".bak-"//ADJUSTL(cp_to_string(j))
END IF
IF (para_env%is_source()) &
INQUIRE (FILE=filename, exist=exist)
CALL para_env%bcast(exist)
IF ((.NOT. exist) .AND. (i < not_read)) THEN
not_read = i
END IF
END DO
IF (not_read == 1) THEN
density_guess = restart_guess
filename = TRIM(file_name)
IF (para_env%is_source()) INQUIRE (FILE=filename, exist=exist)
CALL para_env%bcast(exist)
IF (.NOT. exist) THEN
CALL cp_warn(__LOCATION__, &
"User requested to restart the wavefunction from a series of restart files named: "// &
TRIM(file_name)//" with extensions (.bak-n). These files do not exist."// &
" Even trying to switch to a plain restart wave-function failes because the"// &
" file named: "//TRIM(file_name)//" does not exist. Please check the existence of"// &
" the file or change properly the value of the keyword WFN_RESTART_FILE_NAME."// &
" Calculation continues using ATOMIC GUESS. ")
density_guess = safe_density_guess
END IF
END IF
last_read = not_read - 1
END IF
did_guess = .FALSE.
IF (dft_control%correct_el_density_dip) THEN
tot_corr_zeff = qs_env%total_zeff_corr
!WRITE(*,*) "tot_corr_zeff = ", tot_corr_zeff
IF ((ABS(tot_corr_zeff) > 0.0_dp) .AND. (density_guess /= restart_guess)) THEN
CALL cp_warn(__LOCATION__, &
"Use SCF_GUESS RESTART in conjunction with "// &
"CORE_CORRECTION /= 0.0 and SURFACE_DIPOLE_CORRECTION TRUE. "// &
"It is always advisable to perform SURFACE_DIPOLE_CORRECTION "// &
"after a simulation without the surface dipole correction "// &
"and using the ensuing wavefunction restart file. ")
END IF
END IF
ounit = -1
print_log = .FALSE.
print_history_log = .FALSE.
IF (para_env%is_source()) THEN
CALL section_vals_val_get(dft_section, &
"SCF%PRINT%RESTART%LOG_PRINT_KEY", &
l_val=print_log)
CALL section_vals_val_get(dft_section, &
"SCF%PRINT%RESTART_HISTORY%LOG_PRINT_KEY", &
l_val=print_history_log)
IF (print_log .OR. print_history_log) THEN
ounit = cp_logger_get_default_io_unit(logger)
END IF
END IF
IF (density_guess == restart_guess) THEN
IF (ounit > 0) THEN
WRITE (UNIT=ounit, FMT="(/,T2,A)") &
"WFN_RESTART| Reading restart file"
END IF
IF (do_kpoints) THEN
natoms = SIZE(particle_set)
CALL read_kpoints_restart(rho_ao_kp, kpoints, work1, &
natoms, para_env, id_nr, dft_section, natom_mismatch)
IF (natom_mismatch) density_guess = safe_density_guess
ELSE
CALL read_mo_set_from_restart(mo_array, atomic_kind_set, qs_kind_set, particle_set, para_env, &
id_nr=id_nr, multiplicity=dft_control%multiplicity, dft_section=dft_section, &
natom_mismatch=natom_mismatch, out_unit=ounit)
IF (natom_mismatch) THEN
density_guess = safe_density_guess
ELSE
DO ispin = 1, nspin
IF (scf_control%level_shift /= 0.0_dp) THEN
CALL get_mo_set(mo_set=mo_array(ispin), mo_coeff=mo_coeff)
CALL cp_fm_to_fm(mo_coeff, ortho)
END IF
! make all nmo vectors present orthonormal
CALL get_mo_set(mo_set=mo_array(ispin), &
mo_coeff=mo_coeff, nmo=nmo, homo=homo)
IF (has_unit_metric) THEN
CALL make_basis_simple(mo_coeff, nmo)
ELSEIF (dft_control%smear) THEN
CALL make_basis_lowdin(vmatrix=mo_coeff, ncol=nmo, &
matrix_s=s_sparse(1)%matrix)
ELSE
! ortho so that one can restart for different positions (basis sets?)
CALL make_basis_sm(mo_coeff, homo, s_sparse(1)%matrix)
END IF
! only alpha spin is kept for restricted
IF (dft_control%restricted) EXIT
END DO
IF (dft_control%restricted) CALL mo_set_restrict(mo_array)
IF (.NOT. scf_control%diagonalization%mom) THEN
IF (dft_control%correct_surf_dip) THEN
IF (ABS(tot_corr_zeff) > 0.0_dp) THEN
CALL set_mo_occupation(mo_array, smear=qs_env%scf_control%smear, &
tot_zeff_corr=tot_corr_zeff)
ELSE
CALL set_mo_occupation(mo_array, smear=qs_env%scf_control%smear)
END IF
ELSE
CALL set_mo_occupation(mo_array, smear=qs_env%scf_control%smear)
END IF
END IF
DO ispin = 1, nspin
IF (scf_control%use_ot .OR. scf_env%method == ot_diag_method_nr) THEN !fm->dbcsr
CALL copy_fm_to_dbcsr(mo_array(ispin)%mo_coeff, &
mo_array(ispin)%mo_coeff_b) !fm->dbcsr
END IF !fm->dbcsr
CALL calculate_density_matrix(mo_array(ispin), &
p_rmpv(ispin)%matrix)
END DO
END IF ! natom_mismatch
END IF
! Maximum Overlap Method
IF (scf_control%diagonalization%mom) THEN
CALL do_mom_guess(nspin, mo_array, scf_control, p_rmpv)
END IF
did_guess = .TRUE.
END IF
IF (density_guess == history_guess) THEN
IF (not_read > 1) THEN
IF (ounit > 0) THEN
WRITE (UNIT=ounit, FMT="(/,T2,A)") &
"WFN_RESTART| Reading restart file history"
END IF
DO i = 1, last_read
j = last_read - i
CALL read_mo_set_from_restart(mo_array, atomic_kind_set, qs_kind_set, particle_set, para_env, &
id_nr=j, multiplicity=dft_control%multiplicity, &
dft_section=dft_section, out_unit=ounit)
DO ispin = 1, nspin
IF (scf_control%level_shift /= 0.0_dp) THEN
CALL get_mo_set(mo_set=mo_array(ispin), mo_coeff=mo_coeff)
CALL cp_fm_to_fm(mo_coeff, ortho)
END IF
! make all nmo vectors present orthonormal
CALL get_mo_set(mo_set=mo_array(ispin), mo_coeff=mo_coeff, nmo=nmo, homo=homo)
IF (has_unit_metric) THEN
CALL make_basis_simple(mo_coeff, nmo)
ELSE
! ortho so that one can restart for different positions (basis sets?)
CALL make_basis_sm(mo_coeff, homo, s_sparse(1)%matrix)
END IF
! only alpha spin is kept for restricted
IF (dft_control%restricted) EXIT
END DO
IF (dft_control%restricted) CALL mo_set_restrict(mo_array)
DO ispin = 1, nspin
CALL set_mo_occupation(mo_set=mo_array(ispin), &
smear=qs_env%scf_control%smear)
END DO
DO ispin = 1, nspin
IF (scf_control%use_ot .OR. scf_env%method == ot_diag_method_nr) THEN !fm->dbcsr
CALL copy_fm_to_dbcsr(mo_array(ispin)%mo_coeff, &
mo_array(ispin)%mo_coeff_b) !fm->dbcsr
END IF !fm->dbcsr
CALL calculate_density_matrix(mo_array(ispin), p_rmpv(ispin)%matrix)
END DO
! Write to extrapolation pipeline
CALL wfi_update(wf_history=qs_env%wf_history, qs_env=qs_env, dt=1.0_dp)
END DO
END IF
did_guess = .TRUE.
END IF
IF (density_guess == random_guess) THEN
DO ispin = 1, nspin
CALL get_mo_set(mo_set=mo_array(ispin), &
mo_coeff=mo_coeff, nmo=nmo)
CALL cp_fm_init_random(mo_coeff, nmo)
IF (has_unit_metric) THEN
CALL make_basis_simple(mo_coeff, nmo)
ELSE
CALL make_basis_sm(mo_coeff, nmo, s_sparse(1)%matrix)
END IF
! only alpha spin is kept for restricted
IF (dft_control%restricted) EXIT
END DO
IF (dft_control%restricted) CALL mo_set_restrict(mo_array)
DO ispin = 1, nspin
CALL set_mo_occupation(mo_set=mo_array(ispin), &
smear=qs_env%scf_control%smear)
END DO
DO ispin = 1, nspin
IF (scf_control%use_ot .OR. scf_env%method == ot_diag_method_nr) THEN !fm->dbcsr
CALL copy_fm_to_dbcsr(mo_array(ispin)%mo_coeff, &
mo_array(ispin)%mo_coeff_b) !fm->dbcsr
END IF !fm->dbcsr
CALL calculate_density_matrix(mo_array(ispin), p_rmpv(ispin)%matrix)
END DO
did_guess = .TRUE.
END IF
IF (density_guess == core_guess) THEN
IF (do_kpoints) THEN
CPABORT("calculate_first_density_matrix: core_guess not implemented for k-points")
END IF
owns_ortho = .FALSE.
IF (.NOT. ASSOCIATED(work1)) THEN
need_wm = .TRUE.
CPASSERT(.NOT. ASSOCIATED(work2))
CPASSERT(.NOT. ASSOCIATED(ortho))
ELSE
need_wm = .FALSE.
CPASSERT(ASSOCIATED(work2))
IF (.NOT. ASSOCIATED(ortho)) THEN
ALLOCATE (ortho)
owns_ortho = .TRUE.
END IF
END IF
IF (need_wm) THEN
CALL get_mo_set(mo_set=mo_array(1), mo_coeff=moa)
CALL cp_fm_get_info(moa, matrix_struct=ao_mo_struct)
CALL cp_fm_struct_get(ao_mo_struct, nrow_global=nao, nrow_block=nblocks)
CALL cp_fm_struct_create(fmstruct=ao_ao_struct, &
nrow_block=nblocks, &
ncol_block=nblocks, &
nrow_global=nao, &
ncol_global=nao, &
template_fmstruct=ao_mo_struct)
ALLOCATE (work1(1))
ALLOCATE (work2, ortho)
CALL cp_fm_create(work1(1), ao_ao_struct)
CALL cp_fm_create(work2, ao_ao_struct)
CALL cp_fm_create(ortho, ao_ao_struct)
CALL copy_dbcsr_to_fm(matrix_s_kp(1, 1)%matrix, ortho)
CALL cp_fm_cholesky_decompose(ortho)
CALL cp_fm_struct_release(ao_ao_struct)
END IF
ispin = 1
! Load core Hamiltonian into work matrix
CALL copy_dbcsr_to_fm(h_core_sparse(1)%matrix, work1(ispin))
! Diagonalize the core Hamiltonian matrix and retrieve a first set of
! molecular orbitals (MOs)
IF (has_unit_metric) THEN
CALL eigensolver_simple(matrix_ks=work1(ispin), &
mo_set=mo_array(ispin), &
work=work2, &
do_level_shift=.FALSE., &
level_shift=0.0_dp, &
use_jacobi=.FALSE., jacobi_threshold=0._dp)
ELSE
CALL eigensolver(matrix_ks_fm=work1(ispin), &
mo_set=mo_array(ispin), &
ortho=ortho, &
work=work2, &
cholesky_method=scf_env%cholesky_method, &
do_level_shift=.FALSE., &
level_shift=0.0_dp, &
use_jacobi=.FALSE.)
END IF
! Open shell case: copy alpha MOs to beta MOs
IF (nspin == 2) THEN
CALL get_mo_set(mo_set=mo_array(1), mo_coeff=moa)
CALL get_mo_set(mo_set=mo_array(2), mo_coeff=mob, nmo=nmo)
CALL cp_fm_to_fm(moa, mob, nmo)
END IF
! Build an initial density matrix (for each spin in the case of
! an open shell calculation) from the first MOs set
DO ispin = 1, nspin
CALL set_mo_occupation(mo_set=mo_array(ispin), smear=scf_control%smear)
CALL calculate_density_matrix(mo_array(ispin), p_rmpv(ispin)%matrix)
END DO
! release intermediate matrices
IF (need_wm) THEN
CALL cp_fm_release(ortho)
CALL cp_fm_release(work2)
CALL cp_fm_release(work1(1))
DEALLOCATE (ortho, work2)
DEALLOCATE (work1)
NULLIFY (work1, work2, ortho)
ELSE IF (owns_ortho) THEN
DEALLOCATE (ortho)
END IF
did_guess = .TRUE.
END IF
IF (density_guess == atomic_guess) THEN
subsys_section => section_vals_get_subs_vals(input, "SUBSYS")
ounit = cp_print_key_unit_nr(logger, subsys_section, "PRINT%KINDS", extension=".Log")
IF (ounit > 0) THEN
WRITE (UNIT=ounit, FMT="(/,(T2,A))") &
"Atomic guess: The first density matrix is obtained in terms of atomic orbitals", &
" and electronic configurations assigned to each atomic kind"
END IF
CALL calculate_atomic_block_dm(p_rmpv, s_sparse(1)%matrix, atomic_kind_set, qs_kind_set, &
nspin, nelectron_spin, ounit, para_env)
DO ispin = 1, nspin
! The orbital transformation method (OT) requires not only an
! initial density matrix, but also an initial wavefunction (MO set)
IF (ofgpw .AND. (scf_control%use_ot .OR. scf_env%method == ot_diag_method_nr)) THEN
! get orbitals later
ELSE
IF (need_mos) THEN
IF (dft_control%restricted .AND. (ispin == 2)) THEN
CALL mo_set_restrict(mo_array)
ELSE
CALL get_mo_set(mo_set=mo_array(ispin), &
mo_coeff=mo_coeff, &
nmo=nmo, nao=nao, homo=homo)
CALL cp_fm_set_all(mo_coeff, 0.0_dp)
CALL cp_fm_init_random(mo_coeff, nmo)
CALL cp_fm_create(sv, mo_coeff%matrix_struct, "SV")
! multiply times PS
IF (has_unit_metric) THEN
CALL cp_fm_to_fm(mo_coeff, sv)
ELSE
! PS*C(:,1:nomo)+C(:,nomo+1:nmo) (nomo=NINT(nelectron/maxocc))
CALL cp_dbcsr_sm_fm_multiply(s_sparse(1)%matrix, mo_coeff, sv, nmo)
END IF
CALL cp_dbcsr_sm_fm_multiply(p_rmpv(ispin)%matrix, sv, mo_coeff, homo)
CALL cp_fm_release(sv)
! and ortho the result
IF (has_unit_metric) THEN
CALL make_basis_simple(mo_coeff, nmo)
ELSE
CALL make_basis_sm(mo_coeff, nmo, s_sparse(1)%matrix)
END IF
END IF
CALL set_mo_occupation(mo_set=mo_array(ispin), &
smear=qs_env%scf_control%smear)
CALL copy_fm_to_dbcsr(mo_array(ispin)%mo_coeff, &
mo_array(ispin)%mo_coeff_b) !fm->dbcsr
CALL calculate_density_matrix(mo_array(ispin), &
p_rmpv(ispin)%matrix)
END IF
! adjust el_density in case surface_dipole_correction is switched
! on and CORE_CORRECTION is non-zero
IF (scf_env%method == general_diag_method_nr) THEN
IF (dft_control%correct_surf_dip) THEN
IF (ABS(tot_corr_zeff) > 0.0_dp) THEN
CALL get_mo_set(mo_set=mo_array(ispin), &
mo_coeff=mo_coeff, &
nmo=nmo, nao=nao, homo=homo)
CALL cp_fm_set_all(mo_coeff, 0.0_dp)
CALL cp_fm_init_random(mo_coeff, nmo)
CALL cp_fm_create(sv, mo_coeff%matrix_struct, "SV")
! multiply times PS
IF (has_unit_metric) THEN
CALL cp_fm_to_fm(mo_coeff, sv)
ELSE
! PS*C(:,1:nomo)+C(:,nomo+1:nmo) (nomo=NINT(nelectron/maxocc))
CALL cp_dbcsr_sm_fm_multiply(s_sparse(1)%matrix, mo_coeff, sv, nmo)
END IF
CALL cp_dbcsr_sm_fm_multiply(p_rmpv(ispin)%matrix, sv, mo_coeff, homo)
CALL cp_fm_release(sv)
! and ortho the result
IF (has_unit_metric) THEN
CALL make_basis_simple(mo_coeff, nmo)
ELSE
CALL make_basis_sm(mo_coeff, nmo, s_sparse(1)%matrix)
END IF
CALL set_mo_occupation(mo_set=mo_array(ispin), smear=qs_env%scf_control%smear, &
tot_zeff_corr=tot_corr_zeff)
CALL calculate_density_matrix(mo_array(ispin), &
p_rmpv(ispin)%matrix)
END IF
END IF
END IF
END IF
END DO
IF (ofgpw .AND. (scf_control%use_ot .OR. scf_env%method == ot_diag_method_nr)) THEN
! We fit a function to the square root of the density
CALL qs_rho_update_rho(rho, qs_env)
CPASSERT(1 == 0)
! CALL cp_fm_create(sv,mo_coeff%matrix_struct,"SV")
! DO ispin=1,nspin
! CALL integrate_ppl_rspace(qs%rho%rho_r(ispin),qs_env)
! CALL cp_cfm_solve(overlap,mos)
! CALL get_mo_set(mo_set=mo_array(ispin),&
! mo_coeff=mo_coeff, nmo=nmo, nao=nao)
! CALL cp_fm_init_random(mo_coeff,nmo)
! END DO
! CALL cp_fm_release(sv)
END IF
IF (scf_control%diagonalization%mom) THEN
CALL do_mom_guess(nspin, mo_array, scf_control, p_rmpv)
END IF
CALL cp_print_key_finished_output(ounit, logger, subsys_section, &
"PRINT%KINDS")
did_guess = .TRUE.
END IF
IF (density_guess == sparse_guess) THEN
IF (ofgpw) CPABORT("SPARSE_GUESS not implemented for OFGPW")
IF (.NOT. scf_control%use_ot) CPABORT("OT needed!")
IF (do_kpoints) THEN
CPABORT("calculate_first_density_matrix: sparse_guess not implemented for k-points")
END IF
eps = 1.0E-5_dp
ounit = cp_logger_get_default_io_unit(logger)
natoms = SIZE(particle_set)
ALLOCATE (kind_of(natoms))
ALLOCATE (first_sgf(natoms), last_sgf(natoms))
checksum = dbcsr_checksum(s_sparse(1)%matrix)
i = dbcsr_get_num_blocks(s_sparse(1)%matrix); CALL para_env%sum(i)
IF (ounit > 0) WRITE (ounit, *) 'S nblks', i, ' checksum', checksum
CALL dbcsr_filter(s_sparse(1)%matrix, eps)
checksum = dbcsr_checksum(s_sparse(1)%matrix)
i = dbcsr_get_num_blocks(s_sparse(1)%matrix); CALL para_env%sum(i)
IF (ounit > 0) WRITE (ounit, *) 'S nblks', i, ' checksum', checksum
CALL get_particle_set(particle_set, qs_kind_set, first_sgf=first_sgf, &
last_sgf=last_sgf)
CALL get_atomic_kind_set(atomic_kind_set=atomic_kind_set, kind_of=kind_of)
ALLOCATE (pmat(SIZE(atomic_kind_set)))
rscale = 1._dp
IF (nspin == 2) rscale = 0.5_dp
DO ikind = 1, SIZE(atomic_kind_set)
atomic_kind => atomic_kind_set(ikind)
qs_kind => qs_kind_set(ikind)
NULLIFY (pmat(ikind)%mat)
CALL calculate_atomic_orbitals(atomic_kind, qs_kind, pmat=pmat(ikind)%mat)
NULLIFY (atomic_kind)
END DO
DO ispin = 1, nspin
CALL get_mo_set(mo_set=mo_array(ispin), &
maxocc=maxocc, &
nelectron=nelectron)
!
CALL dbcsr_iterator_start(iter, p_rmpv(ispin)%matrix)
DO WHILE (dbcsr_iterator_blocks_left(iter))
CALL dbcsr_iterator_next_block(iter, irow, icol, pdata, blk)
ikind = kind_of(irow)
IF (icol .EQ. irow) THEN
IF (ispin == 1) THEN
pdata(:, :) = pmat(ikind)%mat(:, :, 1)*rscale + &
pmat(ikind)%mat(:, :, 2)*rscale
ELSE
pdata(:, :) = pmat(ikind)%mat(:, :, 1)*rscale - &
pmat(ikind)%mat(:, :, 2)*rscale
END IF
END IF
END DO
CALL dbcsr_iterator_stop(iter)
!CALL dbcsr_verify_matrix(p_rmpv(ispin)%matrix)
checksum = dbcsr_checksum(p_rmpv(ispin)%matrix)
occ = dbcsr_get_occupation(p_rmpv(ispin)%matrix)
IF (ounit > 0) WRITE (ounit, *) 'P_init occ', occ, ' checksum', checksum
! so far p needs to have the same sparsity as S
!CALL dbcsr_filter(p_rmpv(ispin)%matrix, eps)
!CALL dbcsr_verify_matrix(p_rmpv(ispin)%matrix)
checksum = dbcsr_checksum(p_rmpv(ispin)%matrix)
occ = dbcsr_get_occupation(p_rmpv(ispin)%matrix)
IF (ounit > 0) WRITE (ounit, *) 'P_init occ', occ, ' checksum', checksum
CALL dbcsr_dot(p_rmpv(ispin)%matrix, s_sparse(1)%matrix, trps1)
rscale = REAL(nelectron, dp)/trps1
CALL dbcsr_scale(p_rmpv(ispin)%matrix, rscale)
!CALL dbcsr_verify_matrix(p_rmpv(ispin)%matrix)
checksum = dbcsr_checksum(p_rmpv(ispin)%matrix)
occ = dbcsr_get_occupation(p_rmpv(ispin)%matrix)
IF (ounit > 0) WRITE (ounit, *) 'P occ', occ, ' checksum', checksum
!
! The orbital transformation method (OT) requires not only an
! initial density matrix, but also an initial wavefunction (MO set)
IF (dft_control%restricted .AND. (ispin == 2)) THEN
CALL mo_set_restrict(mo_array)
ELSE
CALL get_mo_set(mo_set=mo_array(ispin), &
mo_coeff=mo_coeff, &
nmo=nmo, nao=nao, homo=homo)
CALL cp_fm_set_all(mo_coeff, 0.0_dp)
n = MAXVAL(last_sgf - first_sgf) + 1
size_atomic_kind_set = SIZE(atomic_kind_set)
ALLOCATE (buff(n, n), sort_kind(size_atomic_kind_set), &
nelec_kind(size_atomic_kind_set))
!
! sort kind vs nbr electron
DO ikind = 1, size_atomic_kind_set
atomic_kind => atomic_kind_set(ikind)
qs_kind => qs_kind_set(ikind)
CALL get_atomic_kind(atomic_kind=atomic_kind, &
natom=natom, &
atom_list=atom_list, &
z=z)
CALL get_qs_kind(qs_kind, nsgf=nsgf, elec_conf=elec_conf, &
basis_set=orb_basis_set, zeff=zeff)
nelec_kind(ikind) = SUM(elec_conf)
END DO
CALL sort(nelec_kind, size_atomic_kind_set, sort_kind)
!
! a -very- naive sparse guess
nmo_tmp = nmo
natoms_tmp = natoms
istart_col = 1
iseed(1) = 4; iseed(2) = 3; iseed(3) = 2; iseed(4) = 1 ! set the seed for dlarnv
DO i = 1, size_atomic_kind_set
ikind = sort_kind(i)
atomic_kind => atomic_kind_set(ikind)
CALL get_atomic_kind(atomic_kind=atomic_kind, &
natom=natom, atom_list=atom_list)
DO iatom = 1, natom
!
atom_a = atom_list(iatom)
istart_row = first_sgf(atom_a)
n_rows = last_sgf(atom_a) - first_sgf(atom_a) + 1
!
! compute the "potential" nbr of states for this atom
n_cols = MAX(INT(REAL(nmo_tmp, dp)/REAL(natoms_tmp, dp)), 1)
IF (n_cols .GT. n_rows) n_cols = n_rows
!
nmo_tmp = nmo_tmp - n_cols
natoms_tmp = natoms_tmp - 1
IF (nmo_tmp .LT. 0 .OR. natoms_tmp .LT. 0) THEN
CPABORT("Wrong1!")
END IF
DO j = 1, n_cols
CALL dlarnv(1, iseed, n_rows, buff(1, j))
END DO
CALL cp_fm_set_submatrix(mo_coeff, buff, istart_row, istart_col, &
n_rows, n_cols)
istart_col = istart_col + n_cols
END DO
END DO
IF (istart_col .LE. nmo) THEN
CPABORT("Wrong2!")
END IF
DEALLOCATE (buff, nelec_kind, sort_kind)
IF (.FALSE.) THEN
ALLOCATE (buff(nao, 1), buff2(nao, 1))
DO i = 1, nmo
CALL cp_fm_get_submatrix(mo_coeff, buff, 1, i, nao, 1)
IF (SUM(buff**2) .LT. 1E-10_dp) THEN
WRITE (*, *) 'wrong', i, SUM(buff**2)
END IF
length = SQRT(DOT_PRODUCT(buff(:, 1), buff(:, 1)))
buff(:, :) = buff(:, :)/length
DO j = i + 1, nmo
CALL cp_fm_get_submatrix(mo_coeff, buff2, 1, j, nao, 1)
length = SQRT(DOT_PRODUCT(buff2(:, 1), buff2(:, 1)))
buff2(:, :) = buff2(:, :)/length
IF (ABS(DOT_PRODUCT(buff(:, 1), buff2(:, 1)) - 1.0_dp) .LT. 1E-10_dp) THEN
WRITE (*, *) 'wrong2', i, j, DOT_PRODUCT(buff(:, 1), buff2(:, 1))
DO ikind = 1, nao
IF (ABS(mo_coeff%local_data(ikind, i)) .GT. 1e-10_dp) THEN
WRITE (*, *) 'c1', ikind, mo_coeff%local_data(ikind, i)
END IF
IF (ABS(mo_coeff%local_data(ikind, j)) .GT. 1e-10_dp) THEN
WRITE (*, *) 'c2', ikind, mo_coeff%local_data(ikind, j)
END IF
END DO
CPABORT("")
END IF
END DO
END DO
DEALLOCATE (buff, buff2)
END IF
!
CALL cp_fm_to_dbcsr_row_template(mo_dbcsr, mo_coeff, s_sparse(1)%matrix)
!CALL dbcsr_verify_matrix(mo_dbcsr)
checksum = dbcsr_checksum(mo_dbcsr)
occ = dbcsr_get_occupation(mo_dbcsr)
IF (ounit > 0) WRITE (ounit, *) 'C occ', occ, ' checksum', checksum
CALL dbcsr_filter(mo_dbcsr, eps)
!CALL dbcsr_verify_matrix(mo_dbcsr)
occ = dbcsr_get_occupation(mo_dbcsr)
checksum = dbcsr_checksum(mo_dbcsr)
IF (ounit > 0) WRITE (ounit, *) 'C occ', occ, ' checksum', checksum
!
! multiply times PS
IF (has_unit_metric) THEN
CPABORT("has_unit_metric will be removed soon")
END IF
!
! S*C
CALL dbcsr_copy(mo_tmp_dbcsr, mo_dbcsr, name="mo_tmp")
CALL dbcsr_multiply("N", "N", 1.0_dp, s_sparse(1)%matrix, mo_dbcsr, &
0.0_dp, mo_tmp_dbcsr, &
retain_sparsity=.TRUE.)
!CALL dbcsr_verify_matrix(mo_tmp_dbcsr)
checksum = dbcsr_checksum(mo_tmp_dbcsr)
occ = dbcsr_get_occupation(mo_tmp_dbcsr)
IF (ounit > 0) WRITE (ounit, *) 'S*C occ', occ, ' checksum', checksum
CALL dbcsr_filter(mo_tmp_dbcsr, eps)
!CALL dbcsr_verify_matrix(mo_tmp_dbcsr)
checksum = dbcsr_checksum(mo_tmp_dbcsr)
occ = dbcsr_get_occupation(mo_tmp_dbcsr)
IF (ounit > 0) WRITE (ounit, *) 'S*C occ', occ, ' checksum', checksum
!
! P*SC
! the destroy is needed for the moment to avoid memory leaks !
! This one is not needed because _destroy takes care of zeroing.
CALL dbcsr_multiply("N", "N", 1.0_dp, p_rmpv(ispin)%matrix, &
mo_tmp_dbcsr, 0.0_dp, mo_dbcsr)
IF (.FALSE.) CALL dbcsr_verify_matrix(mo_dbcsr)
checksum = dbcsr_checksum(mo_dbcsr)
occ = dbcsr_get_occupation(mo_dbcsr)
IF (ounit > 0) WRITE (ounit, *) 'P*SC occ', occ, ' checksum', checksum
CALL dbcsr_filter(mo_dbcsr, eps)
!CALL dbcsr_verify_matrix(mo_dbcsr)
checksum = dbcsr_checksum(mo_dbcsr)
occ = dbcsr_get_occupation(mo_dbcsr)
IF (ounit > 0) WRITE (ounit, *) 'P*SC occ', occ, ' checksum', checksum
!
CALL copy_dbcsr_to_fm(mo_dbcsr, mo_coeff)
CALL dbcsr_release(mo_dbcsr)
CALL dbcsr_release(mo_tmp_dbcsr)