-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathhfbtho_solverV1.f90
executable file
·8188 lines (8171 loc) · 359 KB
/
hfbtho_solverV1.f90
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
!***********************************************************************
!
! Copyright (c) 2016, Lawrence Livermore National Security, LLC.
! Produced at the Lawrence Livermore National
! Laboratory.
! Written by Nicolas Schunck, [email protected]
!
! LLNL-CODE-728299 All rights reserved.
! LLNL-CODE-573953 All rights reserved.
!
! Copyright 2017, R. Navarro Perez, N. Schunck, R. Lasseri, C. Zhang,
! J. Sarich
! Copyright 2012, M.V. Stoitsov, N. Schunck, M. Kortelainen, H.A. Nam,
! N. Michel, J. Sarich, S. Wild
! Copyright 2005, M.V. Stoitsov, J. Dobaczewski, W. Nazarewicz, P.Ring
!
! This file is part of HFBTHO.
!
! HFBTHO is free software: you can redistribute it and/or modify it
! under the terms of the GNU General Public License as published by
! the Free Software Foundation, either version 3 of the License, or
! (at your option) any later version.
!
! HFBTHO is distributed in the hope that it will be useful, but
! WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU General Public License for more details.
!
! You should have received a copy of the GNU General Public License
! along with HFBTHO. If not, see <http://www.gnu.org/licenses/>.
!
! OUR NOTICE AND TERMS AND CONDITIONS OF THE GNU GENERAL PUBLIC
! LICENSE
!
! Our Preamble Notice
!
! A. This notice is required to be provided under our contract
! with the U.S. Department of Energy (DOE). This work was
! produced at the Lawrence Livermore National Laboratory under
! Contract No. DE-AC52-07NA27344 with the DOE.
! B. Neither the United States Government nor Lawrence Livermore
! National Security, LLC nor any of their employees, makes any
! warranty, express or implied, or assumes any liability or
! responsibility for the accuracy, completeness, or usefulness
! of any information, apparatus, product, or process disclosed,
! or represents that its use would not infringe privately-owned
! rights.
! C. Also, reference herein to any specific commercial products,
! process, or services by trade name, trademark, manufacturer
! or otherwise does not necessarily constitute or imply its
! endorsement, recommendation, or favoring by the United States
! Government or Lawrence Livermore National Security, LLC. The
! views and opinions of authors expressed herein do not
! necessarily state or reflect those of the United States
! Government or Lawrence Livermore National Security, LLC, and
! shall not be used for advertising or product endorsement
! purposes.
!
! The precise terms and conditions for copying, distribution and
! modification are contained in the file COPYING.
!
!***********************************************************************
! ==================================================================== !
! !
! DFT SOLVER PACKAGE !
! !
! ==================================================================== !
!----------------------------------------------------------------------
!> This module provides the main \theCode DFT solver. It includes
!> routines for the calculation and diagonalization of the HFB matrix;
!> the definition of densities on the quadrature mesh and in configuration
!> space; the self-consistent loop the calculation of expectation values of
!> observables. It also includes a package of a few routines to perform
!> particle number projection.
!>
!> @author
!> Mario Stoitsov, Nicolas Schunck, Markus kortelainen, Rodrigo Navarro Perez
!----------------------------------------------------------------------
! Subroutines: - HFBTHO_DFT_SOLVER
! - heading
! - thodefh
! - preparer
! - coordinateLST
! - iter
! - hfbdiag
! - ALambda
! - Canonical
! - resu
! - initialize_HFBTHO_NAMELIST
! - read_HFBTHO_NAMELIST
! - check_consistency
! - initialize_HFBTHO_SOLVER
! - base0
! - base
! - gaupol
! - start
! - printRHO
! - gfv
! - sdiag
! - nucleus
! - stab
! - ord
! - tracesln
! - tracesln_qp
! - densitln
! - coulom1
! - coulom
! - coulom_test
! - HartreeDir
! - optHFBTHO
! - DENSIT
! - field
! - gamdel
! - gamdel_gogny
! - expectpj
! - densitpj
! - coulompj
! - broyden_min
! - expect
! - Constraint_or_not
! - getLagrange
! - requested_blocked_level
!----------------------------------------------------------------------
Module HFBTHO_solver
Use HFBTHO_utilities
Use HFBTHO
Use HFBTHO_multipole_moments
Use HFBTHO_collective
Use HFBTHO_fission_fragments
Use HFBTHO_gogny
Use HFBTHO_io
Implicit None
Contains
!-----------------------------------------------------------------------------------------------
!
! Axially-deformed configurational constrained and/or unconstrained Hartree-Fock-Bogoliubov
! calculations with Skyrme-like functionals and delta pairing using the Harmonic-Oscillator
! (HO), and/or Transformed HO (THO) basis with or without reflection symmetry imposed, with
! or without the Lipkin-Nogami procedure. The solver can handle all Skyrme-like functionals,
! DME-functionals, Fayans-functionals, calculate infinite nuclear matter properties, finite
! nuclei (even-even, odd-even, odd-odd), and neutron drops, isoscalar and isovector monopo-
! le FAM QRPA calculations for spherical and deformed nuclei.
!
! All necessary input variables contain the suffix _INI. Below, the complete list of these
! variables with some example values:
!
! ======== hfbtho_NAMELIST.dat
! n00_INI=20; npr1_INI=70; npr2_INI=50; kindhfb_INI=-1; inin_INI=-1
! b0_INI=2.234776; q_INI=0.0; cdef_INI=0.0; cqad_INI=0.5; skyrme_INI='SLY4'; nkblo_INI=0
! ILST_INI=0; keypj_INI=1; iproj_INI=0; npr1pj_INI=0;
! icou_INI=2; IDEBUG_INI=0; npr2pj_INI=0;
! Parity_INI=.False.; epsi_INI=0.00001_pr; MAX_ITER_INI=101
! Add_Pairing_INI=.False.; DO_FITT_INI=.False.; Print_PTHO_Namelist_INI=.True.
!
! ======== from read_UNEDF_NAMELIST
! DMEORDER=-1; DMELDA=0; use_TMR_pairing=0
! HBZERO=20.73553000000000; E2CHARG=1.4399784085965135; CRHO(0)=-933.3423749999999; CRHO(1)=830.0524855000001;
! CDRHO(0)=861.0625000000000; CDRHO(1)=-1064.2732500000; CTAU(0)=57.12868750000000; CTAU(1)=24.65673650000000;
! CRDR(0)=-76.99620312499999; CRDR(1)=15.65713512500000; CRDJ(0)=-92.25000000000000; CRDJ(1)=-30.7500000000000;
! CJ(0)=17.20961150000000; CJ(1)=64.57581250000000; CPV0(0)=-258.2000000000000; CPV0(1)=-258.2000000000000;
! CPV1(0)=0.5000000000000000; CPV1(1)=0.500000000000000; SIGMA=0.1666666666666667; CEXPAR=1.000000000000000;
! E_NM=-15.97214914144462; K_NM=229.9009644826037; SMASS_NM =1.439546988976078; RHO_NM =0.1595387567117334;
! ASS_NM =32.00430281505202; LASS_NM=45.96175148046161; VMASS_NM =1.249838547196253;
! MPI=0.6995945261023822; GA=1.290000000000000; FPI=0.4683223517486062; C1=-0.1598130000000000;
! C3 =-0.6708200000000; C4 =0.6708200000000000; CD =-2.062000000000000; CE=-0.6250000000000;
! LAMBDAX =3.547896604156107; USE_INM=.false.; USE_CM_COR =.true.; USE_DME3N_TERMS=.true.;
! USE_J2TERMS =.true.; USE_CHARGE_DENSITY=.false.; PRINT_NAMELIST=.true.;
!
! Memo:
! - inin_INI switches scratch unconstrained (inin=1,2,3) or constrained (inin=100,200,300)
! calculations. Unconstrained mode begins with a small number of constrained iterations.
! - inin_INI switches unconstrained (inin=-1,-2,-3) or constrained (inin=-100,-200,-300)
! calculations from a previous solution if the latter exists. If not, the solver sets
! inin=Abs(inin) and resumes from scratch.
! - The same holds for odd nuclei. If even-even solution for the odd nucleus does not exists
! it is calculated first.
! - Print_Screen=T/F for n00_INI=+/-: output is generated and written in thoout.dat file
! only if n00_INI>0. For n00_INI<0, the number of shells is set to abs(n00_INI) but all
! output is supressed.
!
! - At the end of the solution, the solver provides all results in the arrays
!
! nucname,ereslbl(1:2),eres(1:ierest)
!
! which contain:
!
! LBL,BLKN,BLKZ,Jsi,JININ,A,N',Z,Efn,Efp,JEtot,Jbett,Jbetn,Jbetp,JQt,JQn,JQp
! JpEn,JpEp,JpDn,JpDp,JAsn,JAsp,Jrt,Jrn,Jrp,Jrc,Jht,Jhn,Jhp,Jqht,Jqhn,Jqhp,
! JKINt,JKINn,JKINp,JSO,JCDIR,JCEX,JDisn,JDisp,JV2Mn,JV2Mp,JILST,JKIND,JL,
! JECMPAV1,JECMPAV2,JECMPAV3,JA,JN,JZ,ITER,UEtot,Ubett,Ubetn,Ubetp,UQt,UQn,
! UQp,Uln,Ulp,UpEn,UpEp,UpDn,UpDp,UAsn,UAsp,Urt,Urn,Urp,Urc,Uht,Uhn,Uhp,
! Uqht,Uqhn,Uqhp,UKINT,UKINN,UKINP,USO,UCDIR,UCEX,UDisn,UDisp,UV2Mn,UV2Mp,
! UECMT,UECMN,UECMP,UROTT,UROTN,UROTP,USQUJT,USQUJN,USQUJP,UCRANT,UCRANN,
! UCRANP,UERIGT,UERIGN,UERIGP,EHFBLN,EHFB,LNbet,LNben,LNbep,LNQt,LNQn,LNQp,
! LNpEn,LNpEp,LNpDn,LNpDp,LNrt,LNrn,LNrp,LNrC,LNam2n,LNam2p,LNe2n,LNe2p,
! BlEqpN,BlDEqpN,BlOvrN,BlEqpZ,BlDEqpZ,BlOvrZ
!
! - Standard inputs from 'hfbtho_NAMELIST.dat'
! - USer-defined functional read from 'UNEDF_NAMELIST.DAT'
! - Final solution stored to files '*.hel' and/or '*.tel'
! - Output is written to files 'thoout.dat', 'thodef.dat' and 'hodef.dat'
! - Output files *.dat may exist as 'thoout','thores','hodenp',
! 'thodenp','thodef','thoene',
! 'thoprc','dat0.1.2.3.4'
! - External accuracy pr/ipr always set in UNEDF module
!
! - n00 Number of oscillator shells
! n00>0 prints to thoout.dat & screen
! n00<0 no print at all
! n00=0 program stops (NB!)
! - b0 Oscillator Basis parameter b0>0 (If b0<0 it takes a default value)
! - beta0 Value of Basis deformation parameter
! - AN Number of neutrons N
! - AZ Number of protons Z
! - FTST' Fayance forces label
! - kind Kind of calculations 1: noLN, -1:LN
! - inin Unconstraint Iterations from peviouse solution
! -1: (from spherical *.hel or *.tel file)
! -2: (from prolate *.hel or *.tel file)
! -3: (from oblate *.hel or *.tel file)
! Unconstraint Iterations fom scratch with
! preliminary constraint at deformation Cbeta, (i):
! 1: Spherical scratch
! 2: Prolate scratch
! 3: Oblate scratch
! Constrained calculation (icstr) at Cbeta, see (i)
! 100, 200, 300 fom scratch
! -100,-200,-300 fom previouse solution
! - blNeutrons: a group responsible for blocking a particular neutron level
! The group consists of 5 numbers, e.g., for 7-[ 3, 0, 3]: 7 -1 3 0 3
! k1 2 \times \Omega
! =0: the whole group (k) is disregarded (n0 blocking)
! >0: blocking in N+1 nucleus
! <0: blocking in N-1 nucleus
! k2 parity (+1 or -1); NB! when k2=0, the ground state walker is applied
! k3,k4,k5 Nilson quantum numbers
! - blProtons: exactly the same as (k) but for protons
!
!-----------------------------------------------------------------------------------------------
Subroutine HFBTHO_DFT_SOLVER
Use HFBTHO_THO
Implicit None
Integer(ipr) :: iw,ib,j,i,it,l,maxi0,icstr0,iterMax,icons,il,kickoff,iexit,ncons_eff
Real(pr) :: epsi0,qq,f,f1,f2,f3,r,g,g1
!-------------------------------------------------------------
! Initializing all according to *_INI values
!-------------------------------------------------------------
Call initialize_HFBTHO_SOLVER
If(ierror_flag.Ne.0) Return
#if(DO_MASSTABLE==1 || DRIP_LINES==1 || DO_PES==1)
If(lout.le.lfile) Open(lfile,file='thoout'//row_string//'.dat',status='unknown')
#else
#ifdef SUPRESS_OUTPUT
If(lout.le.lfile) Open(lfile,file='thoout.dat',status='unknown')
#else
If(lout.lt.lfile) Open(lfile,file='thoout.dat',status='unknown')
#endif
#endif
Call Constraint_or_not(inin_INI,inin,icstr)
If(ierror_flag.Ne.0) Return
!-------------------------------------------------------------------------
! Loop recalculating eventually the even-even solution for an odd nucleus
!-------------------------------------------------------------------------
irestart=0; iexit=0
matrix_elements_calculated = .false.
Do
n00=Abs(n00_INI); b0=b0_INI; q=q_INI; iLST=iLST_INI;
maxi=MAX_ITER_INI; npr(1)=npr_INI(1); npr(2)=npr_INI(2); npr(3)=npr(1)+npr(2);
skyrme=skyrme_INI; kindhfb=kindhfb_INI
keypj=keypj_INI; iproj=iproj_INI; npr1pj=npr1pj_INI; npr2pj=npr2pj_INI;
nkblo=nkblo_INI;
basis_HFODD=basis_HFODD_INI
!-------------------------------------------------------------
! Define the set of constraints
!-------------------------------------------------------------
numberCons=0; kickoff=0
Do l=1,lambdaMax
If(lambda_active(l).Gt.0) numberCons = numberCons + 1
If(lambda_active(l).Lt.0) kickoff = kickoff + 1
End Do
! Add constraint on the neck
If(neck_constraints) Then
numberCons = numberCons + 1
neckLag = zero
!kickoff = kickoff + 1
End If
!
ncons_eff = numberCons + kickoff
If(.Not.Allocated(multLag)) Allocate(multLag(1:lambdaMax)); multLag=zero
If(.Not.Allocated(multLambda)) Allocate(multLambda(1:ncons_eff)); multLambda=0
If(.Not.Allocated(multRequested)) Allocate(multRequested(0:lambdaMax)); multRequested=zero
!
icons=0
Do l=1,lambdaMax
If(lambda_active(l).Gt.0) Then
icons=icons+1
multLambda(icons)=lambda_values(l)
End If
multRequested(l) = expectation_values(l)
End Do
If(neck_constraints) Then
icons=icons+1
multLambda(icons)=0
End If
!-------------------------------------------------------------
! Blocking
!-------------------------------------------------------------
Do it=1,2
If(nkblo(it,1).Ne.0) Then
If(nkblo(it,1).Gt.0) Then
! particle state
npr(it)=npr(it)+1
iparenti(it)=-1
Else
! hole state
npr(it)=npr(it)-1
iparenti(it)=+1
End If
nkblo(it,1)=Abs(nkblo(it,1))
If(nkblo(it,2).Eq.0) Then
! ground state walker
keyblo(it)=keyblo(it)+1 !nkblo(it,1)
If(keyblo(it).Eq.blomax(it)) irestart=0
Else
irestart=0
End If
End If
End Do
!-------------------------------------------------------------
! HFB+HO calculations
!-------------------------------------------------------------
If(ILST.Le.0) Then
icacou=0; icahartree=0
Call preparer(.True.)
If(ierror_flag.Ne.0) Return
! Reading input file: if not possible, start from scratch
Call inout(1,iexit)
If(iexit>0) Then
Call start
Else
Call gamdel(.false.,.false.)
End If
If(ierror_flag.Ne.0) Return
!-------------------------------------------------------------
! Preliminary constrained calculations
!-------------------------------------------------------------
If(kickoff.Gt.0) Then
icstr0=icstr; epsi0=epsi; ! remember accuracy
icstr=1 ! constraint true
epsi=1.0_pr ! small accuracy
if(is_nedf) then
iterMax = maxi; maxi = 25
else
iterMax = maxi; maxi = 10
endif
numberCons=0
Do l=1,lambdaMax
If(Abs(lambda_active(l)).Gt.0) Then
numberCons=numberCons+1
multLambda(numberCons)=lambda_values(l)
End If
End Do
Do iw=lout,lfile
If(Parity) Then
Write(iw,'(/,a,i3,a,i2,a,/)') ' ### INITIAL STAGE(constrained calculations, reflection symmetry used)'
Else
Write(iw,'(/,a,i3,a,i2,a,/)') ' ### INITIAL STAGE(constrained calculations, no reflection symmetry used)'
End If
End Do
Call iter(.True.) ! small constraint iterations
If(ierror_flag.Ne.0) Return
! For the next phase, use only true constraints
icstr=icstr0; epsi=epsi0
maxi = iterMax
numberCons=0
Do l=1,lambdaMax
If(lambda_active(l).Gt.0) Then
numberCons=numberCons+1
multLambda(numberCons)=lambda_values(l)
End If
End Do
End If
!-------------------------------------------------------------
! REGULAR HFB+HO ITERATIONS
!-------------------------------------------------------------
Do iw=lout,lfile
If(Parity) Then
Write(iw,'(/,a,i3,a,i2,a,/)') ' ### REGULAR STAGE (reflection symmetry imposed)'
Else
Write(iw,'(/,a,i3,a,i2,a,/)') ' ### REGULAR STAGE (no reflection symmetry imposed)'
End If
End Do
Call iter(.True.)
If(ierror_flag.Ne.0) Return
Call resu(1)
If(ierror_flag.Ne.0) Return
End If
!! write LST function on disk
!Open(unit=66,file='LST.dat',status='unknown')
!Write(66,'("#",15X,"R",20X,"f(R)",20X,"f^(1)",20X,"f^(2)",20X,"f^(3)")')
!Do il=1,170
! qq=Real(il-1)/10.0_pr*1.0_pr
! If(il.Eq.1) Call thofun(0,g,f,f1,f2,f3,g1,.True.,.True.)
! Call thofun(1,qq,f,f1,f2,f3,r,.False.,.True.)
! Write(66,'(6E24.10)') r,qq,f1,f2,f3
!End Do
!Close(66)
!-------------------------------------------------------------
! HFB+THO calculations from HFB+HO
!-------------------------------------------------------------
If(ILST.Lt.0) Then
ILST1=1; icacou=0; icahartree=0
Call coordinateLST(.False.) ! THO basis
If(ierror_flag.Ne.0) Return
Call densit ! THO densities
If(ierror_flag.Ne.0) Return
Call field ! Nuclear fields
If(ierror_flag.Ne.0) Return
Call iter(.True.) ! HFB+THO iterations
If(ierror_flag.Ne.0) Return
Call resu(1) ! print/record results
If(ierror_flag.Ne.0) Return
End If
!-------------------------------------------------------------
! HFB+THO calculations from *.tel
!-------------------------------------------------------------
If(ILST.Gt.0) Then
If(inin.Gt.0) Then
ierror_flag=ierror_flag+1
ierror_info(ierror_flag)= ' Forbidden iLST>0, inin>0 '
Return
End If
icacou=0; icahartree=0
Call preparer(.True.)
If(ierror_flag.Ne.0) Return
! Reading input file: if not possible, start from scratch
Call inout(1,iexit)
If(iexit>0) Then
Call start
Else
Call gamdel(.false.,.false.)
End If
If(ierror_flag.Ne.0) Return
Call iter(.True.) ! HFB+THO iterations
If(ierror_flag.Ne.0) Return
Call resu(1) ! print/record results
If(ierror_flag.Ne.0) Return
End If
!-------------------------------------------------------------
! Go for the requested blocking state in a case of odd nuclei
! if restarted due to corrupted/missing previous solution
!-------------------------------------------------------------
inin=-Abs(inin)
If(irestart.Eq.0) Exit
End Do
!
End Subroutine HFBTHO_DFT_SOLVER
!=======================================================================
!> Print heading to screen 'lout' and to tape thoout.dat 'lfile'
!=======================================================================
Subroutine heading
#if(USE_OPENMP==1)
Use omp_lib
#endif
#if(USE_MPI==1)
use mpi
#endif
Implicit None
#if(USE_MPI==1)
Integer(ipr) :: rank,ierr
#endif
Integer(ipr) :: iw,idt(8),numThreads,idThread
Character(len=12) rcl(3)
Character(len=50) today
!
#if(USE_OPENMP==1)
!$OMP PARALLEL SHARED(iw,lout,lfile) PRIVATE(numThreads,idThread)
numThreads = omp_get_num_threads()
idThread = omp_get_thread_num()
#if(DO_MASSTABLE==1 || DO_PES==1)
If (idThread .Eq. 0.and.irow.le.1) Then
#elif(DRIP_LINES==1)
If (idThread .Eq. 0.and.irow.le.1.and.calc_counter.eq.1) Then
#elif(USE_MPI==1)
call MPI_COMM_RANK(MPI_COMM_WORLD, rank, ierr)
If (idThread .Eq. 0 .and. rank .eq. 0) Then
#else
If (idThread .Eq. 0) Then
#endif
Do iw=lout,lfile
Write(iw,'("Multi-threading framework with OpenMP:",i2," threads/task")') numThreads
End Do
End If
!$OMP END PARALLEL
#endif
Call Date_and_time(rcl(1),rcl(2),rcl(3),idt)
Write(today,'(a,i2,a,i2,a,i4,a,i2,a,i2,a)')'(',idt(2),'/',idt(3),'/',idt(1),', ',idt(5),':',idt(6),')'
Do iw=lout,lfile
Write(iw,*)
Write(iw,'(" ==========================================")')
Write(iw,'(" FORTRAN 95 CODE (KIND=",i2,") ")') pr
Write(iw,'(" Version: ",a)') Version
Write(iw,'(" ==========================================")')
Write(iw,'(" AXIALLY DEFORMED CONFIGURATIONAL ")')
Write(iw,'(" HARTREE-FOCK-BOGOLIUBOV CALCULATIONS ")')
Write(iw,'(" WITH ")')
Write(iw,'(" SKYRME+DELTA PAIRING OR GOGNY EDFs ")')
Write(iw,'(" USING THE ")')
Write(iw,'(" HARMONIC-OSCILLATOR ")')
Write(iw,'(" AND/OR ")')
Write(iw,'(" TRANSFORMED HARMONIC-OSCILLATOR ")')
Write(iw,'(" BASIS ")')
Write(iw,'(" --- ")')
Write(iw,'(" v1.66 (2005) ")')
Write(iw,'(" Stoitsov, Dobaczewski, Nazarewicz, Ring ")')
Write(iw,'(" v2.00d (2012) ")')
Write(iw,'(" Stoitsov, Schunck, Kortelainen ")')
Write(iw,'(" v3.00 (2016) ")')
Write(iw,'(" Schunck, Navarro Perez ")')
Write(iw,'(" ==========================================")')
Write(iw,'(" Nucleus: ",a," (A=",i4,", N=",i3,", Z=",i3,")")') nucname,npr(1)+npr(2),npr(1),npr(2)
If(Parity) Then
Write(iw,'(" Reflection Symmetry Imposed ")')
Else
Write(iw,'(" No Reflection Symmetry Imposed ")')
End If
Write(iw,'(" ",a)') today
Write(iw,'(" ==========================================")')
Write(iw,*)
Write(iw,*)
End Do
End Subroutine heading
!=======================================================================
!> Print labels to hodef.dat or/and thodef.dat files
!=======================================================================
Subroutine thodefh(iw1)
Implicit None
Integer(ipr) :: iw1
hlabels(1)='LBL'; hlabels(11)='JEtot'; hlabels(21)='JpDp';
hlabels(2)='BLKN'; hlabels(12)='Jbett'; hlabels(22)='JAsn';
hlabels(3)='BLKZ'; hlabels(13)='Jbetn'; hlabels(23)='JAsp';
hlabels(4)='Jsi'; hlabels(14)='Jbetp'; hlabels(24)='Jrt';
hlabels(5)='JININ'; hlabels(15)='JQt'; hlabels(25)='Jrn';
hlabels(6)='A'; hlabels(16)='JQn'; hlabels(26)='Jrp';
hlabels(7)='N'; hlabels(17)='JQp'; hlabels(27)='Jrc';
hlabels(8)='Z'; hlabels(18)='JpEn'; hlabels(28)='Jht';
hlabels(9)='Efn'; hlabels(19)='JpEp'; hlabels(29)='Jhn';
hlabels(10)='Efp'; hlabels(20)='JpDn'; hlabels(30)='Jhp';
!
hlabels(31)='Jqht'; hlabels(41)='JDisp'; hlabels(51)='JN';
hlabels(32)='Jqhn'; hlabels(42)='JV2Mn'; hlabels(52)='JZ';
hlabels(33)='Jqhp'; hlabels(43)='JV2Mp'; hlabels(53)='ITER';
hlabels(34)='JKINt'; hlabels(44)='JILST'; hlabels(54)='UEtot';
hlabels(35)='JKINn'; hlabels(45)='JKIND'; hlabels(55)='Ubett';
hlabels(36)='JKINp'; hlabels(46)='JL'; hlabels(56)='Ubetn';
hlabels(37)='JSO'; hlabels(47)='JECMPAV1'; hlabels(57)='Ubetp';
hlabels(38)='JCDIR'; hlabels(48)='JECMPAV2'; hlabels(58)='UQt';
hlabels(39)='JCEX'; hlabels(49)='JECMPAV3'; hlabels(59)='UQn';
hlabels(40)='JDisn'; hlabels(50)='JA'; hlabels(60)='UQp';
!
hlabels(61)='Uln'; hlabels(71)='Urp'; hlabels(81)='UKINP';
hlabels(62)='Ulp'; hlabels(72)='Urc'; hlabels(82)='USO';
hlabels(63)='UpEn'; hlabels(73)='Uht'; hlabels(83)='UCDIR';
hlabels(64)='UpEp'; hlabels(74)='Uhn'; hlabels(84)='UCEX';
hlabels(65)='UpDn'; hlabels(75)='Uhp'; hlabels(85)='UDisn';
hlabels(66)='UpDp'; hlabels(76)='Uqht'; hlabels(86)='UDisp';
hlabels(67)='UAsn'; hlabels(77)='Uqhn'; hlabels(87)='UV2Mn';
hlabels(68)='UAsp'; hlabels(78)='Uqhp'; hlabels(88)='UV2Mp';
hlabels(69)='Urt'; hlabels(79)='UKINT'; hlabels(89)='UECMT';
hlabels(70)='Urn'; hlabels(80)='UKINN'; hlabels(90)='UECMN';
!
hlabels(91)='UECMP'; hlabels(101)='UERIGT'; hlabels(111)='LNQp';
hlabels(92)='UROTT'; hlabels(102)='UERIGN'; hlabels(112)='LNpEn';
hlabels(93)='UROTN'; hlabels(103)='UERIGP'; hlabels(113)='LNpEp';
hlabels(94)='UROTP'; hlabels(104)='EHFBLN'; hlabels(114)='LNpDn';
hlabels(95)='USQUJT'; hlabels(105)='EHFB'; hlabels(115)='LNpDp';
hlabels(96)='USQUJN'; hlabels(106)='LNbet'; hlabels(116)='LNrt';
hlabels(97)='USQUJP'; hlabels(107)='LNben'; hlabels(117)='LNrn';
hlabels(98)='UCRANT'; hlabels(108)='LNbep'; hlabels(118)='LNrp';
hlabels(99)='UCRANN'; hlabels(109)='LNQt'; hlabels(119)='LNrC';
hlabels(100)='UCRANP'; hlabels(110)='LNQn'; hlabels(120)='LNam2n';
!
hlabels(121)='LNam2p';
hlabels(122)='LNe2n';
hlabels(123)='LNe2p';
hlabels(124)='BlEqpN';
hlabels(125)='BlDEqpN';
hlabels(126)='BlOvrN';
hlabels(127)='BlEqpZ';
hlabels(128)='BlDEqpZ';
hlabels(129)='BlOvrZ';
!
Write(iw1,'((1x,a,2x),6x,660(a,2x))') hlabels
!
! HELP
!Do i=1,129
! Write(iw1,'(1x,i3,a,a)',advance='NO') i,':',trim(hlabels(i))
!End Do
! 1:LBL 2:BLKN 3:BLKZ 4:Jsi 5:JININ 6:A 7:N 8:Z 9:Efn 10:Efp
! 11:JEtot 12:Jbett 13:Jbetn 14:Jbetp 15:JQt 16:JQn 17:JQp 18:JpEn 19:JpEp 20:JpDn
! 21:JpDp 22:JAsn 23:JAsp 24:Jrt 25:Jrn 26:Jrp 27:Jrc 28:Jht 29:Jhn 30:Jhp
! 31:Jqht 32:Jqhn 33:Jqhp 34:JKINt 35:JKINn 36:JKINp 37:JSO 38:JCDIR 39:JCEX 40:JDisn
! 41:JDisp 42:JV2Mn 43:JV2Mp 44:JILST 45:JKIND 46:JL 47:JECMPAV1 48:JECMPAV2 49:JECMPAV3 50:JA
! 51:JN 52:JZ 53:ITER 54:UEtot 55:Ubett 56:Ubetn 57:Ubetp 58:UQt 59:UQn 60:UQp
! 61:Uln 62:Ulp 63:UpEn 64:UpEp 65:UpDn 66:UpDp 67:UAsn 68:UAsp 69:Urt 70:Urn
! 71:Urp 72:Urc 73:Uht 74:Uhn 75:Uhp 76:Uqht 77:Uqhn 78:Uqhp 79:UKINT 80:UKINN
! 81:UKINP 82:USO 83:UCDIR 84:UCEX 85:UDisn 86:UDisp 87:UV2Mn 88:UV2Mp 89:UECMT 90:UECMN
! 91:UECMP 92:UROTT 93:UROTN 94:UROTP 95:USQUJT 96:USQUJN 97:USQUJP 98:UCRANT 99:UCRANN 100:UCRANP
! 101:UERIGT 102:UERIGN 103:UERIGP 104:EHFBLN 105:EHFB 106:LNbet 107:LNben 108:LNbep 109:LNQt 110:LNQn
! 111:LNQp 112:LNpEn 113:LNpEp 114:LNpDn 115:LNpDp 116:LNrt 117:LNrn 118:LNrp 119:LNrC 120:LNam2n
! 121:LNam2p 122:LNe2n 123:LNe2p 124:BlEqpN 125:BlDEqpN 126:BlOvrN 127:BlEqpZ 128:BlDEqpZ 129:BlOvrZ
End Subroutine thodefh
!=======================================================================
!> Allocates arrays depending on the nunber of shells/states and on
!> the number of Gauss points
!=======================================================================
Subroutine thoalloc
Implicit None
Integer :: ier,ib,ND
!
! number of int.points
If(Parity) Then
ngh=ngh_INI; ngl=ngl_INI; nleg=nleg_INI !Yesp
Else
ngh=2*ngh_INI; ngl=ngl_INI; nleg=nleg_INI !Nop
End If
!
!nbx=2*n00+1 ! maximal number of k-blocks
!ntx=(n00+1)*(n00+2)*(n00+3)/6 ! max.num. p/n levels
!nzx=n00 ! maximal nz-quantum number
!nrx=n00/2+1 ! maximal nr-quantum number
!nlx=n00 ! maximal ml-quantum number
!ndx=(n00+2)*(n00+2)/4 ! maximal dim. of one k-block
!nhhdim=number of nonzero HH matrix elements
!
nzrlx=(nzx+1)*(nrx+1)*(nlx+1) ! phy(:,:,nzrlx)
nghl=ngh*ngl ! nghl=ngh*ngl
nqx=ndx*ndx; nb2x=nbx+nbx; ndx2=ndx+ndx
ilnqx=ilpj*nqx; ilnghl=ilpj*nghl
nhfbx=ndx+ndx; nhfbqx=nhfbx*nhfbx; nkx=ntx; ndxs=ndx*(ndx+1)/2
!-----------------------------------------
!Arrays depending on gauss points
!-----------------------------------------
If(Allocated(xleg)) Deallocate(xleg,wleg)
If(nleg.Gt.0) Allocate(xleg(nleg),wleg(nleg))
If(Allocated(xh)) Deallocate(xh,wh,xl,sxl,wl,vc ,vhbn,vn,vrn,vzn,vdn,vsn,dvn,vhbp,vp,vrp,vzp,vdp,vsp,dvp, &
vSZFIn,vSFIZn,vSRFIn,vSFIRn,vSZFIp,vSFIZp,vSRFIp,vSFIRp, &
fl,fli,fh,fd,fp1,fp2,fp3,fp4,fp5,fp6, fs1,fs2,fs3,fs4,fs5,fs6, &
wdcor,wdcori,cou,vDHartree,vhart00,vhart01,vhart11)
Allocate(xh(ngh),wh(ngh),xl(ngl),sxl(ngl),wl(ngl),vc(nghl,nghl))
Allocate(vhbn(nghl),vn(nghl),vrn(nghl),vzn(nghl),vdn(nghl),vsn(nghl),dvn(nghl), &
vhbp(nghl),vp(nghl),vrp(nghl),vzp(nghl),vdp(nghl),vsp(nghl),dvp(nghl), &
vSZFIn(nghl),vSFIZn(nghl),vSRFIn(nghl),vSFIRn(nghl), &
vSZFIp(nghl),vSFIZp(nghl),vSRFIp(nghl),vSFIRp(nghl))
Allocate(fl(nghl),fli(nghl),fh(nghl),fd(nghl),fp1(nghl),fp2(nghl),fp3(nghl), &
fp4(nghl),fp5(nghl),fp6(nghl),fs1(nghl),fs2(nghl),fs3(nghl),fs4(nghl), &
fs5(nghl),fs6(nghl),wdcor(nghl),wdcori(nghl),cou(nghl),vDHartree(nghl,2), &
vhart00(nghl,nghl),vhart01(nghl,nghl),vhart11(nghl,nghl))
If(Allocated(aka)) Deallocate(aka,ro,tau,dro,dj,NABLAR,NABLAZ,SZFI,SFIZ,SRFI,SFIR)
Allocate(aka(nghl,2),ro(nghl,2),tau(nghl,2),dro(nghl,2),dj(nghl,2), &
SZFI(nghl,2),SFIZ(nghl,2),SRFI(nghl,2),SFIR(nghl,2),NABLAR(nghl,2),NABLAZ(nghl,2))
If(Allocated(qfield)) Deallocate(qfield)
Allocate(qfield(nghl,lambdaMax+1)) ! constraining fields: lambdaMax multipoles + neck
If(Allocated(MEFFn)) Deallocate(MEFFn,MEFFp)
Allocate(MEFFn(nghl),MEFFp(nghl))
If(Allocated(geff_inv)) Deallocate(geff_inv)
Allocate(geff_inv(nghl,2))
!-----------------------------------------
! Arrays depending on configurations
!-----------------------------------------
If(Allocated(rk)) Deallocate(rk,ak,qh,qh1,ql,ql1,nz,nr,nl,ns,npar,id &
,ia,ikb,ipb,ka,kd,tb,txb,numax,ek,dk,vk,vk1,uk,vkmax,ddc,ddc1,hfb1,lcanon)
Allocate(rk(nqx,nb2x),ak(nqx,nb2x),qh(0:nzx,1:ngh+1) &
,qh1(0:nzx,1:ngh+1),ql(0:nrx,0:nlx,1:ngl+1),ql1(0:nrx,0:nlx,1:ngl+1) &
,nz(ntx),nr(ntx),nl(ntx),ns(ntx),npar(ntx),id(nbx),ia(nbx),ikb(nbx),lcanon(0:nbx,2) &
,ipb(nbx),ka(nbx,2),kd(nbx,2),tb(ntx),txb(nbx),numax(0:nkx,2) &
,ek(nkx,2),dk(nkx,2),vk(nkx,2),vk1(nkx,2),uk(nkx,2),vkmax(nkx,2) &
,ddc(ndx,nkx,2),ddc1(ndx,nkx,2),hfb1(nhfbx,2))
If(finite_range.or.coulomb_gaussian) Then
If(Allocated(nrr)) Deallocate(nrr,nll,nss,noo,nzz,nzzx)
Allocate(nrr(ntx),nll(ntx),nss(ntx),noo(ntx),nzz(ntx,ntx),nzzx(ntx))
End If
!-----------------------------------------
! HFB Arrays
!-----------------------------------------
If(Allocated(erhfb)) Deallocate(erhfb,drhfb,erhfb1,drhfb1)
Allocate(erhfb(nkx),drhfb(nkx),erhfb1(nkx),drhfb1(nkx))
If(Allocated(hfb)) Deallocate(hfb,zhfb,evvk,hfbcan,evvkcan)
Allocate(hfb(ndx2,ndx2),zhfb(ndx2),evvk(ndx2),hfbcan(ndx,ndx),evvkcan(ndx))
If(Allocated(AN)) Deallocate(AN,ANk,PFIU,PFID,FIU,FID,FIUR,FIDR,FIUD2N,FIDD2N,FIUZ,FIDZ)
Allocate(AN(nqx),ANk(nqx),PFIU(ndx),PFID(ndx),FIU(ndx),FID(ndx) &
,FIUR(ndx),FIDR(ndx),FIUD2N(ndx),FIDD2N(ndx),FIUZ(ndx),FIDZ(ndx))
!-----------------------------------------
! Optimal LAPACK storage
!-----------------------------------------
If(Allocated(alwork)) Deallocate(alwork,lwork)
#if(SWITCH_ESSL==0)
ialwork=1+6*ndx+2*ndx**2; ilwork=3+5*ndx;
Allocate(alwork(ialwork),lwork(ilwork));alwork = 0.0; lwork = 1
#else
ialwork=0; ilwork=5*ndx;
Allocate(alwork(1),lwork(ilwork));alwork = 0.0; lwork = 0
#endif
!ialwork=1; ilwork=1;
!If(Allocated(alwork)) Deallocate(alwork,lwork)
!Allocate(alwork(ialwork),lwork(ilwork))
!ier=0; Call DSYEVD('V','L',ndx2,hfb,ndx2,evvk,ALWORK,-1,LWORK,-1,ier)
!If(ier.Ne.0) Then
! ierror_flag=ierror_flag+1
! ierror_info(ierror_flag)='STOP: FATAL ERROR CONDITION IN DSYEVD'
! Return
!End If
!ialwork=Int(alwork(1)); ilwork=lwork(1)
!If(Allocated(alwork)) Deallocate(alwork,lwork)
!Allocate(alwork(ialwork),lwork(ilwork))
!-----------------------------------------
! Eqp, U,V
!-----------------------------------------
If(Allocated(RVqpN)) Deallocate(RVqpN,RVqpP,RUqpN,RUqpP,REqpN,REqpP)
Allocate(RVqpN(nuv),RVqpP(nuv),RUqpN(nuv),RUqpP(nuv),REqpN(nqp),REqpP(nqp))
If(Allocated(KpwiP)) Deallocate(KpwiP,KpwiN,KqpN,KqpP)
Allocate(KpwiN(nqp),KpwiP(nqp),KqpN(nqp),KqpP(nqp))
If(Allocated(fn_T)) Deallocate(fn_T,fp_T)
Allocate(fn_T(nqp),fp_T(nqp))
!-----------------------------------------
! PNP ARRAYS: CONF. AND GAUGE ANGLE
!-----------------------------------------
If(Allocated(exp1iphy))Deallocate(ropj,taupj,dropj,djpj,akapj,coupj,pjk &
,SZFIpj,SFIZpj,SRFIpj,SFIRpj,epj,cpj,ypj,rpj,ddepj,phypj,sinphy &
,exp1iphy,exp2iphy,exp1iphym,exp2iphym)
Allocate(ropj(nghl,ilpj,2),taupj(nghl,ilpj,2),dropj(nghl,ilpj,2) &
,djpj(nghl,ilpj,2),akapj(nghl,ilpj,2),coupj(nghl,ilpj),pjk(ilpj,2) &
,SZFIpj(nghl,ilpj,2),SFIZpj(nghl,ilpj,2),SRFIpj(nghl,ilpj,2) &
,SFIRpj(nghl,ilpj,2),epj(ilpj,2),cpj(nkx,ilpj,2),ypj(nkx,ilpj,2) &
,rpj(nkx,ilpj,2),ddepj(nqx,ilpj,nb2x),phypj(ilpj),sinphy(ilpj), &
exp1iphy(ilpj),exp2iphy(ilpj),exp1iphym(ilpj),exp2iphym(ilpj))
!-----------------------------------------
! FIELDS INITIALIZATION (NB! optimize)
!-----------------------------------------
ro=zero; tau=zero; dro=zero; dj=zero; aka=zero; rk=zero
vn=zero; vsn=zero; vhbn=zero; vrn=zero; vzn=zero; vdn=zero;
vp=zero; vsp=zero; vhbp=zero; vrp=zero; vzp=zero; vdp=zero;
dvn=zero; dvp=zero;
vSFIZn=zero; vSZFIn=zero; vSFIRn=zero; vSRFIn=zero; vDHartree=zero;
vSFIZp=zero; vSZFIp=zero; vSFIRp=zero; vSRFIp=zero;
! Jason
If(Allocated(allhfb)) Then
Do ib=1,oldnb
Deallocate(allhfb(ib)%arr,allevvk(ib)%arr,allalwork(ib)%arr,alllwork(ib)%arr)
End Do
Deallocate (allhfb,allevvk,allalwork,alllwork)
Deallocate (allIALWORK,allILWORK,allISUPPZ)
End If
If (Allocated(allibro)) Deallocate(allibro)
!
End Subroutine thoalloc
!=======================================================================
!> Setup routine: initializes most variables - including NAMELISTS; set
!> up the basis and the quadrature grid; computes matrix elements of the
!> Gogny force (if finite range); prints summary information of the run
!=======================================================================
Subroutine preparer(lpr)
Use HFBTHO_gauss
Implicit None
Logical :: lpr
Integer(ipr) :: iw,l,icount
Real(pr) :: amas_base
!
If(n00.Eq.0) Then
ierror_flag=ierror_flag+1
ierror_info(ierror_flag)=' STOP: No more nuclei pass to the solver'
Return
End If
!-----------------------------------------
! select the symbol of the nucleus
!-----------------------------------------
Call nucleus(1,npr(2),nucname)
If(ierror_flag.Ne.0) Return
!-----------------------------------------
! print headings to screen/'thoout.dat'
!-----------------------------------------
If(lpr) Then
Call heading
Call print_functional_parameters()
Do iw=lout,lfile
If(ierror_flag.Ne.0) Return
If(Print_HFBTHO_Namelist) Then
Write(iw,'(100(2x,a,f15.8))')
Write(iw,'(100(2x,a,f15.8))') 'NAMELIST CONTENT (copy/past to hfbtho_NAMELIST.dat and modify)'
Write(iw,'(100(2x,a,f15.8))') '-------------------------------------------------------------'
Write(iw,HFBTHO_GENERAL)
Write(iw,HFBTHO_INITIAL)
Write(iw,HFBTHO_ITERATIONS)
Write(iw,HFBTHO_FUNCTIONAL)
Write(iw,HFBTHO_PAIRING)
Write(iw,HFBTHO_CONSTRAINTS)
Write(iw,HFBTHO_BLOCKING)
Write(iw,HFBTHO_PROJECTION)
Write(iw,HFBTHO_FEATURES)
Write(iw,HFBTHO_TDDFT)
Write(iw,HFBTHO_NECK)
Write(iw,HFBTHO_TEMPERATURE)
Write(iw,HFBTHO_DEBUG)
End If
End Do
End If
!-----------------------------------------
! particle number as real variable
!-----------------------------------------
tz(1)=Real(npr(1),Kind=pr); tz(2)=Real(npr(2),Kind=pr)
If(fragment_properties) Then
tz=tz_fragments(1:2)
End If
amas=tz(1)+tz(2)
! In case of blocking, force the same oscillator length for the core
! even-even nucleus
amas_base = amas
If(.Not.odd_noBlock) Then
If(nkblo_INI(1,1).Gt.0.And.npr(1).Eq.2*(npr(1)/2)) amas_base = amas_base+one
If(nkblo_INI(1,1).lt.0.And.npr(1).Eq.2*(npr(1)/2)) amas_base = amas_base-one
End If
If(.Not.odd_noBlock) Then
If(nkblo_INI(2,1).Gt.0.And.npr(2).Eq.2*(npr(2)/2)) amas_base = amas_base+one
If(nkblo_INI(2,1).lt.0.And.npr(2).Eq.2*(npr(2)/2)) amas_base = amas_base-one
End If
drhoi=zero
!-----------------------------------------
! default combinations
!-----------------------------------------
chargee2=e2charg
coex=-chargee2*(three/pi)**p13; cex=-0.750_pr*coex
!--------------------------------------------------------
! hbzero from forces [hqc**2/(two*amu)] EOedit for SV-min
!--------------------------------------------------------
if(.not.hb0_charge_dependent) then
hbzeron = hbzero
hbzerop = hbzero
! hb0=hbzero; If (use_cm_cor) hb0=hb0*(one-one/amas)
! hb0n=hbzeron; If (use_cm_cor) hb0n=hb0n*(one-one/amas)
! hb0p=hbzerop; If (use_cm_cor) hb0p=hb0p*(one-one/amas)
!else if () then
!STUFF
else
! hbzero = 0.5_pr*(hbzeron + hbzerop)
hb0=hbzero; If (use_cm_cor) hb0=hb0*(one-one/amas)
hb0n=hbzeron; If (use_cm_cor) hb0n=hb0n*(one-one/amas)
hb0p=hbzerop; If (use_cm_cor) hb0p=hb0p*(one-one/amas)
endif
hb0=hbzero; If (use_cm_cor) hb0=hb0*(one-one/amas)
hb0n=hbzeron; If (use_cm_cor) hb0n=hb0n*(one-one/amas)
hb0p=hbzerop; If (use_cm_cor) hb0p=hb0p*(one-one/amas)
!-----------------------------------------
! basis parameter q
!-----------------------------------------
beta0=q; q=Exp((3.0_pr*Sqrt(5.0_pr/(16.0_pr*pi)))*beta0)
!-----------------------------------------
! basis parameters b0,bp,bz
!-----------------------------------------
If(b0.Le.zero) Then
! define oscillator frequency from default with empirical factor 1.2,
! and set length accordingly
r00=r0*amas_base**p13; r02=r00**2; r04=r02**2
hom=41.0_pr*amas_base**(-p13)*r0
b0=Sqrt(two*hbzero/hom)
Else
! define oscillator frequency from user-defined length, and set default
! empirical factor accordingly
hom=hqc**2/(amn*b0**2)
r0=(hom/41.0_pr)*amas**(p13)
r00=r0*amas**p13; r02=r00**2; r04=r02**2
End If
bp=b0*q**(-one/6.0_pr); bz=b0*q**(one/3.0_pr); bpp=bp*bp
!-----------------------------------------
! constraint in terms of beta
!-----------------------------------------
ty20=Sqrt(5.0_pr/pi)*hom/b0**2/two
!-----------------------------------------
! projection: number of grid points
!-----------------------------------------
keypj=Max(1,keypj); ilpj=keypj; ilpj2=ilpj**2;
If(iproj.Eq.0) Then
npr1pj=npr(1); npr2pj=npr(2)
Else
npr1pj=npr(1)+npr1pj; npr2pj=npr(2)+npr2pj
End If
!-----------------------------------------
! blocking window
!-----------------------------------------
!pwiblo=Min(Max(25.0_pr/Sqrt(Real(npr(1)+npr(2),Kind=pr)),2.0_pr),8.0_pr)
pwiblo=1.0_pr
!-----------------------------------------
! THO
!-----------------------------------------
ass=zero; iasswrong=0
!-----------------------------------------
! iterations
!-----------------------------------------
etot=zero; varmas=zero; rms=zero; ept=-two; del=one; alast=-seven; siold=one
varmasNZ=zero; pjmassNZ=zero; ass=zero; skass=zero
!---------------------------------------------------------
! statistics to screen('lout')/file('lfile')
!---------------------------------------------------------
If(lpr) Then
Do iw=lout,lfile
Write(iw,*)
Write(iw,'(a)') ' ---------------------------------------'
Write(iw,'(a)') ' Characteristics of the run '
Write(iw,'(a)') ' ---------------------------------------'
Write(iw,'(a,i5)') ' Output file ................: ',lfile
Write(iw,'(a,2x,a2,i4)') ' Nucleus ....................: ',nucname,npr(1)+npr(2)
Write(iw,'(a,i5)') ' Number of HO shells ........: ',n00
Write(iw,'(a,f20.14)') ' HO length b0 (fm) ..........: ',b0
Write(iw,'(a,f8.3,a,f8.3)') ' Basis deformation ..........: beta0=',beta0,' q=',q
Write(iw,'(a,5(1x,e15.8))') ' HO: b0,1/b0,bp,bz,q ........: ',b0,one/b0,bp,bz,q
Write(iw,'(a,3(1x,e15.8))') ' h**2/(2m_n), cmc, e**2 .....: ',hbzeron,hb0n,chargee2
Write(iw,'(a,3(1x,e15.8))') ' h**2/(2m_p), cmc, e**2 .....: ',hbzerop,hb0p,chargee2
Write(iw,'(a,2(1X,e15.8))') ' hom=f*41.0_pr*A^{-1/3}, f...: ',hom,r0
If(automatic_basis) Then
Write(iw,'(a)') ' Adjusted basis is ........: ON'
End If
If(iLST.Eq.0) Then ! HFB+HO case only
iLST1=0
Write(iw,'(a)') ' THO basis is ...............: OFF'
Else ! HFB+THO case
Write(iw,'(a)') ' THO basis is ...............: ON'
If(iLST.Gt.0) Then ! HFB+THO only
iLST1=1
If(inin.Gt.0) Then
ierror_flag=ierror_flag+1
ierror_info(ierror_flag)=' Stop: Forbidden iLST>0, inin>0 combination.'
Return
End If
Write(iw,'(a)') ' THO parameters from tholst.wel'
Else ! HFB+THO after HFB+HO
iLST1=0
Write(iw,'(a)') ' HFB+THO after a HFB+HO run '
End If
End If
Write(iw,'(a,i5)') ' Maximal number of iterations: ',maxi
Write(iw,'(a,f6.3)') ' Initial mixing parameter ...: ',xmix
If(inin.Eq.1) Then
Write(iw,'(a)') ' Initial w.f. ...............: from spherical scratch'
End If
If(inin.Eq.2) Then
Write(iw,'(a)') ' Initial w.f. ...............: from prolate scratch'
End If
If(inin.Eq.3) Then
Write(iw,'(a)') ' Initial w.f. ...............: from oblate scratch'
End If
If(inin.Lt.0) Then
Write(iw,'(a)') ' Initial wave functions from : tape'
End If
Write(iw,'(a,3x,a)') ' Energy functional ..........: ',skyrme
If(finite_range) Then
Write(iw,'(a,f6.2,a)') ' with a finite-range central force'
End If
If(icou.Eq.-4) Write(iw,'(a)') ' direct Coulomb by sum of Gaussians, exchange Coulomb with Slater approximation'
If(icou.Eq.-3) Write(iw,'(a)') ' direct, exchange and pairing Coulomb by sum of Gaussians'
If(icou.Eq.-2) Write(iw,'(a)') ' direct and exchange Coulomb by sum of Gaussians'
If(icou.Eq.-1) Write(iw,'(a)') ' direct Coulomb force only by sum of Gaussians'
If(icou.Eq. 0) Write(iw,'(a)') ' without Coulomb forces'
If(icou.Eq. 1) Write(iw,'(a)') ' direct Coulomb by substitution method'
If(icou.Eq. 2) Write(iw,'(a)') ' direct Coulomb by substitution method, exchange Coulomb with Slater approximation'
If(kindhfb.Lt.0) Then
Write(iw,'(a)') ' Lipkin-Nogami procedure is .: ON'
Else
Write(iw,'(a)') ' Lipkin-Nogami procedure is .: OFF'
End If
If(ilpj-1.Eq.0) Then
Write(iw,'(a)') ' PAV procedure is ...........: OFF'
Else
Write(iw,'(a)') ' PAV procedure is ...........: ON'
Write(iw,'(a,i5)') ' Number of gauge points....: ',keypj
End If
If(icstr.Eq.0) Then
Write(iw,'(a)') ' Constraint calculation is ..: OFF'
Else
Write(iw,'(a)') ' Constraint calculation is ..: ON'
icount=0
Do l=1,8
If(Abs(lambda_active(l)).Gt.0) Then
icount=icount+1
Write(iw,'(a,i1,a,i1,a,f8.3)') ' Constraint ',icount,' .............: lambda=',l, &
' Ql=',multRequested(l)
End If
End Do
If(neck_constraints) Then
icount=icount+1
Write(iw,'(a,i1,a,a,f8.3)') ' Neck ',icount,' .............: lambda=0', &
' Ql=',neckRequested
End If
End If
If(keyblo(1).Ne.0) Then