-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathm_paw_gaussfit.F90
2215 lines (1957 loc) · 55.2 KB
/
m_paw_gaussfit.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
!!****m* ABINIT/m_paw_gaussfit
!! NAME
!! m_paw_gaussfit
!!
!! FUNCTION
!! Module to fit PAW related data to sums of gaussians
!!
!! COPYRIGHT
!! Copyright (C) 2012-2022 ABINIT group (T. Rangel)
!! This file is distributed under the terms of the
!! GNU General Public License, see ~abinit/COPYING
!! or http://www.gnu.org/copyleft/gpl.txt .
!!
!! NOTES
!! FOR DEVELOPPERS: in order to preserve the portability of libPAW library,
!! please consult ~abinit/src/??_libpaw/libpaw-coding-rules.txt
!!
!! SOURCE
#include "libpaw.h"
module m_paw_gaussfit
USE_DEFS
USE_MSG_HANDLING
USE_MPI_WRAPPERS
USE_MEMORY_PROFILING
use m_paw_numeric, only : paw_splint, paw_spline
use m_pawrad, only : pawrad_type, pawrad_init, pawrad_deducer0, pawrad_free, pawrad_ifromr
implicit none
private
public:: gaussfit_projector !fit non-local projectors to a sum of gaussians
public:: gaussfit_main !main routine to fit
!Routines related to MPI:
private:: gaussfit_mpi_set_weight
private:: gaussfit_mpi_remove_item
private:: gaussfit_mpi_add_item
private:: gaussfit_mpi_assign
private:: gaussfit_mpi_main
private:: gaussfit_mpi_calc_deviation
private:: gaussfit_mpi_swap
!Routines related to fitting:
private:: gaussfit_fit !fit a function to gaussians
private:: gaussfit_calc_deriv_r !calc. derivatives for real gauss.
private:: gaussfit_calc_deriv_c !calc. derivatives for cplex. gauss. of form 1
private:: gaussfit_calc_deriv_c2 !calc. derivatives for cplex. gauss. of form 2
private:: gaussfit_calc_deriv_c3 !calc. derivatives for cplex. gauss. of form 3
private:: gaussfit_calc_deriv_c4 !calc. derivatives for cplex. gauss. of form 4
private:: gaussfit_rlsf !retreined least squares fit
private:: gaussfit_chisq_alpha_beta
!set parameters for LSF (for 5 different forms of gauss. sums):
private:: gaussfit_set_param1
private:: gaussfit_set_param2
private:: gaussfit_set_param3
private:: gaussfit_set_param4
private:: gaussfit_set_param5
private:: gaussfit_constrains_init !initialize constrains
private:: gaussfit_apply_constrains !apply cons. to get new params.
!!***
integer,private,parameter:: positive=2
integer,private,parameter:: restricted=3
integer,private,parameter:: restricted_and_positive=4
CONTAINS
!===========================================================
!!***
!!****f* m_paw_gaussfit/gaussfit_main
!! NAME
!! gaussfit_main
!!
!! FUNCTION
!! Fits a given input function f(r) to a sum of gaussians
!!
!! INPUTS
!! nterm_bounds= sets the minimum and maximum number of terms to be taken into account.
!! mparam= maximum number of parameters
!! if(option==1)mparam=nterm_bounds(2)*4
!! if(option==2)mparam=nterm_bounds(2)*6
!! if(option==3)mparam=nterm_bounds(2)*2
!! if(option==4)mparam=nterm_bounds(2)*4
!! nr= number of real space points
!! pawrad= pawrad type
!! option=1 fit to a1 cos(a2 x^2)+ a3 sin( a4 x^2)
!! 2 fit to a1 exp(-a2 x^2)*(a3 cos (a4 x^2) + a5 sin (a6 x^2) )
!! 3 fit to a1 cos (k x^2) + a2 sin (k x^2)
!! 4 fit to a1 exp(-a2 x^2)* (a3 cos(k x^2)+ a4 sin (k x^2))
!! Given these definitions, the number of complex gaussians are:
!! ngauss=4*nterm (for option=1,2)
!! ngauss=2*nterm (for option=3,4)
!! outfile= filename to write out fitted functions.
!! rpaw=paw radius
!! y(nr)= function to fit
!!
!! OUTPUT
!! nparam_out= number of parameters found.
!! param_out(nparam_out)= parameters (coefficients and factors of complex gaussians).
!!
!! SOURCE
subroutine gaussfit_main(mparam,nparam_out,nterm_bounds,nr,&
& param_out,pawrad,option,outfile,rpaw,y,comm_mpi)
!Arguments ------------------------------------
integer,intent(in)::mparam,nr,nterm_bounds(2)
integer,intent(in)::option
integer,intent(in),optional:: comm_mpi
real(dp),intent(in)::rpaw
real(dp),intent(inout)::y(nr)
character(80),intent(in)::outfile
type(pawrad_type),intent(in) :: pawrad
integer,intent(out)::nparam_out
real(dp),intent(out)::param_out(mparam)
!Local variables-------------------------------
!scalars
logical,parameter::modify_y=.false. !used only for plotting purposes
integer :: ichisq,ierr,ii,jj
integer :: master,me
integer :: maxiter,minterm,my_chisq_size,counts_all
integer :: ngauss,nparam,nproc,nterm
integer :: verbosity
real(dp) :: chisq,chisq_min
!real(dp) :: T1,T2 !uncomment for timming
!arrays
integer::constrains(mparam)
integer::proc_dist(nterm_bounds(1):nterm_bounds(2))
integer,allocatable::counts(:),disp(:),map_nterm(:)
real(dp)::limit(mparam),param_tmp(mparam),weight(mparam)
real(dp),allocatable::chisq_array(:),recv_buf(:),send_buf(:)
real(dp),allocatable::y_out(:)
character(len=500) :: msg
! *************************************************************************
!initialize variables
maxiter=200
!
!initialize mpi quantities:
master=0; me=0; nproc=1; proc_dist=1;
if(present(comm_mpi)) then
me=xpaw_mpi_comm_rank(comm_mpi)
nproc=xpaw_mpi_comm_size(comm_mpi)
end if
if(nproc>1) then
! Find distribution (master)
if(me==master) then
call gaussfit_mpi_main(nproc,nterm_bounds,proc_dist)
end if
! send distribution to all processors
call xpaw_mpi_bcast(proc_dist(nterm_bounds(1):nterm_bounds(2)),&
& master,comm_mpi,ierr)
end if
!Set size of chisq treated by each proc
my_chisq_size=nterm_bounds(2)-nterm_bounds(1)+1 !all terms
if(nproc>1 .and. .not. me==master) then
do nterm=nterm_bounds(1),nterm_bounds(2)
if(.not. proc_dist(nterm)==me+1) cycle
my_chisq_size=my_chisq_size+1
end do
end if
!
!Allocate objects
!
LIBPAW_ALLOCATE(y_out,(nr))
LIBPAW_ALLOCATE(chisq_array,(my_chisq_size))
if(master==me ) then
LIBPAW_BOUND1_ALLOCATE(map_nterm,BOUNDS(nterm_bounds(1),nterm_bounds(2)))
jj=1
do ii=1,nproc
do nterm=nterm_bounds(1),nterm_bounds(2)
if(proc_dist(nterm)==ii) then
map_nterm(nterm)=jj
jj=jj+1
end if
end do
end do
end if
!
!fill with zeros
!
chisq_array=zero
nparam_out=0
y_out=0.d0
verbosity=1 !print the minimum at the terminal
!
ichisq=0
do nterm=nterm_bounds(1),nterm_bounds(2)
! mpi distribution
if(.not. proc_dist(nterm)==me+1) cycle
ichisq=ichisq+1
! call CPU_TIME(T1)
if(option==1) nparam=nterm*4
if(option==2) nparam=nterm*6
if(option==3) nparam=nterm*2
if(option==4) nparam=nterm*4
! set initial guess
! if(option==1) then
! call gaussfit_set_param3(nterm,nparam,param_tmp(1:nparam),sep(ii))
! elseif(option==2) then
! call gaussfit_set_param1(nterm,nparam,nr,&
!& param_tmp(1:nparam),sep(ii),pawrad%rad(1:nr),y)
if(option==3) then
call gaussfit_set_param4(nparam,param_tmp(1:nparam))
elseif(option==4) then
call gaussfit_set_param5(nterm,nparam,nr,&
& param_tmp(1:nparam),rpaw,y)
end if
!
!
call gaussfit_constrains_init(weight(1:nparam),constrains(1:nparam),&
& limit(1:nparam),nparam,nterm,nr,option,rpaw,y)
!
call gaussfit_fit(chisq,constrains(1:nparam),&
& limit(1:nparam),maxiter,nparam,nterm,nr,option,outfile,param_tmp,&
& verbosity,weight(1:nparam),pawrad%rad(1:nr),y,y_out)
! if there was an error, set chisq to very high
! if there is an NaN, for instance:
if(abs(chisq+1.d0)<tol8) chisq=99999
if(abs(chisq)==chisq*chisq) chisq=99999
if(chisq .ne. chisq) chisq=99999
!
chisq_array(ichisq)=chisq
! call CPU_TIME(T2)
! print *, 'Time for fit ', T2-T1, 'seconds.'
end do !nterm
!mpicast results
!send distribution to master
if(nproc>1) then
! Prepare communications:
LIBPAW_ALLOCATE(counts,(nproc))
counts(:)=0
do nterm=nterm_bounds(1),nterm_bounds(2)
counts(proc_dist(nterm))=counts(proc_dist(nterm))+1
end do
counts_all=sum(counts)
LIBPAW_ALLOCATE(send_buf,(counts(me+1)))
send_buf(:)=chisq_array(1:counts(me+1))
if(me==master) then
LIBPAW_ALLOCATE(recv_buf,(counts_all))
else
LIBPAW_ALLOCATE(recv_buf,(1))
end if
LIBPAW_ALLOCATE(disp,(nproc))
disp(1)=0
do ii=2,nproc
disp(ii)=disp(ii-1)+counts(ii-1)
end do
! communicate all info to master
call xpaw_mpi_gatherv(send_buf,counts(me+1),recv_buf,&
& counts,disp,master,comm_mpi,ierr)
! fill in chisq_array with all received info:
if(master==me) then
do ii=1,counts_all
chisq_array(ii)=recv_buf(ii)
end do
end if
! Deallocate MPI arrays:
LIBPAW_DEALLOCATE(recv_buf)
LIBPAW_DEALLOCATE(counts)
LIBPAW_DEALLOCATE(disp)
LIBPAW_DEALLOCATE(send_buf)
end if
!Print out info:
if(me==master) then
write(msg,'(3a)')'Preliminary results (with only 200 iter.):',ch10,' ngauss chisq'
call wrtout(std_out,msg,'COLL')
do nterm=nterm_bounds(1),nterm_bounds(2)
if(option==1) ngauss=nterm*4
if(option==2) ngauss=nterm*4
if(option==3) ngauss=nterm*2
if(option==4) ngauss=nterm*2
write(msg,'(i4,2x,e13.6,1x)')ngauss,&
& chisq_array(map_nterm(nterm))
call wrtout(std_out,msg,'COLL')
end do
end if
!get minterm for best accuracy:
if(me==master) then
chisq_min=999999999
do nterm=nterm_bounds(1),nterm_bounds(2)
if(chisq_array(map_nterm(nterm))<chisq_min) then
chisq_min=chisq_array(map_nterm(nterm))
minterm=nterm
end if
end do
!run again with the best parameters
nterm=minterm
if(option==1)nparam=4*nterm
if(option==2)nparam=6*nterm
if(option==3)nparam=2*nterm
if(option==4)nparam=4*nterm
! set initial guess
if (option==3) then
call gaussfit_set_param4(nparam,param_tmp(1:nparam))
elseif (option==4) then
call gaussfit_set_param5(nterm,nparam,nr,&
& param_tmp(1:nparam),rpaw,y)
end if
!
call gaussfit_constrains_init(weight(1:nparam),constrains(1:nparam),&
& limit(1:nparam),nparam,nterm,nr,option,rpaw,y)
!
verbosity=1;
maxiter=1000 !this time we do more iterations
call gaussfit_fit(chisq,constrains(1:nparam),&
& limit(1:nparam),maxiter,nparam,nterm,nr,option,outfile,param_tmp,&
& verbosity,weight(1:nparam),pawrad%rad(1:nr),y,y_out)
! Write out best solution
write(msg,'(3a)')"Best solution (with more iterations):",ch10," ngauss chisq"
call wrtout(std_out,msg,'COLL')
if(option==1) ngauss=nterm*4
if(option==2) ngauss=nterm*4
if(option==3) ngauss=nterm*2
if(option==4) ngauss=nterm*2
write(msg,'(i4,2x,e13.6,1x)')ngauss,chisq
call wrtout(std_out,msg,'COLL')
end if
!Fill output variables
!and communicate results to all procs:
if(option==1) then
nparam_out=minterm*4
elseif (option==2) then
nparam_out=minterm*6
elseif (option==3) then
nparam_out=minterm*2
elseif (option==4) then
nparam_out=minterm*4
end if
!communicate
if(nproc>1) then
call xpaw_mpi_bcast(nparam_out,master,comm_mpi,ierr)
call xpaw_mpi_bcast(param_tmp(1:nparam_out),master,comm_mpi,ierr)
end if !nproc>1
param_out(:)=param_tmp
if(modify_y) then
! call xpaw_mpi_scatterv(y_out,nr,mpi_displ,y_out,nr,0,comm_mpi,ierr)
y=y_out !at output modify y for the fitted y
end if
LIBPAW_DEALLOCATE(y_out)
LIBPAW_DEALLOCATE(chisq_array)
if(me==master) then
LIBPAW_DEALLOCATE(map_nterm)
end if
end subroutine gaussfit_main
!!***
!----------------------------------------------------------------------
!!****f* m_paw_gaussfit/gaussfit_mpi_set_weight
!! NAME
!! gaussfit_mpi_set_weight
!!
!! FUNCTION
!! It sets a weight to the number of
!! Gaussians used.
!! This was calculated by measuring the time
!! it takes to fit a projector with different
!! number of gaussians.
!!
!! INPUTS
!!
!! OUTPUT
!!
!! SOURCE
subroutine gaussfit_mpi_set_weight(f,x)
!Arguments ------------------------------------
integer,intent(in)::x
integer,intent(out)::f
!Local variables ------------------------------
real(dp)::a,b,c,d,ff,xx
!************************************************************************
!The following parameters were obtained
!from the time (in seconds)
!it takes to fit a given projector
!using from 1 to 100 gaussians.
a = -0.374137d0 ! +/- 2.013 (538.1%)
b = 0.207854d0 ! +/- 0.3385 (162.8%)
c = 0.0266371d0 ! +/- 0.01534 (57.59%)
d = 0.000152476d0 ! +/- 0.0001978 (129.7%)
xx=real(x,dp)
ff=a+b*xx+c*xx**2+d*xx**3
f=max(1,ceiling(ff))
end subroutine gaussfit_mpi_set_weight
!!***
!----------------------------------------------------------------------
!!****f* m_paw_gaussfit/gaussfit_mpi_remove_item
!! NAME
!! gaussfit_mpi_remove_item
!!
!! FUNCTION
!!
!! INPUTS
!!
!! OUTPUT
!!
!! SOURCE
subroutine gaussfit_mpi_remove_item(iterm,pload)
!Arguments ------------------------------------
integer,intent(in)::iterm
integer,intent(inout)::pload
!Local variables ------------------------------
integer:: f_i
!***********************************************************************
call gaussfit_mpi_set_weight(f_i,iterm)
pload=pload-f_i
end subroutine gaussfit_mpi_remove_item
!!***
!----------------------------------------------------------------------
!!****f* m_paw_gaussfit/gaussfit_mpi_add_item
!! NAME
!! gaussfit_mpi_add_item
!!
!! FUNCTION
!!
!! INPUTS
!!
!! OUTPUT
!!
!! SOURCE
subroutine gaussfit_mpi_add_item(iterm,pload)
!Arguments ------------------------------------
integer,intent(in)::iterm
integer,intent(inout)::pload
!Local variables ------------------------------
integer:: f_i
!************************************************************************
call gaussfit_mpi_set_weight(f_i,iterm)
pload=pload+f_i
end subroutine gaussfit_mpi_add_item
!!***
!----------------------------------------------------------------------
!!****f* m_paw_gaussfit/gaussfit_mpi_calc_deviation
!! NAME
!! gaussfit_mpi_calc_deviation
!!
!! FUNCTION
!!
!! INPUTS
!!
!! OUTPUT
!!
!! SOURCE
subroutine gaussfit_mpi_calc_deviation(deviation,nproc,proc_load)
!Arguments ------------------------------------
integer,intent(in)::nproc
integer,intent(in)::proc_load(nproc)
integer,intent(out)::deviation
!Local variables ------------------------------
integer:: jproc,kproc,kload,jload
!************************************************************************
! deviation=0
! do jproc=1,nproc
! deviation=deviation+abs(proc_load(jproc)-ideal)
! end do
deviation=0
do jproc=1,nproc
jload=proc_load(jproc)
do kproc=1,jproc
kload=proc_load(kproc)
deviation=deviation+abs(kload-jload)
end do
end do
end subroutine gaussfit_mpi_calc_deviation
!!***
!----------------------------------------------------------------------
!!****f* m_paw_gaussfit/gaussfit_mpi_swap
!! NAME
!! gaussfit_mpi_swap
!!
!! FUNCTION
!!
!! INPUTS
!!
!! OUTPUT
!!
!! SOURCE
subroutine gaussfit_mpi_swap(iterm,jterm,&
& nproc,nterm_bounds,proc_dist,proc_load)
!Arguments ------------------------------------
integer,intent(in)::iterm,jterm,nproc,nterm_bounds(2)
integer,intent(inout)::proc_dist(nterm_bounds(1):nterm_bounds(2)),proc_load(nproc)
!Local variables ------------------------------
integer:: deviation1,deviation2
integer:: iproc,jproc
!************************************************************************
!Calculate initial state
call gaussfit_mpi_calc_deviation(deviation1,nproc,proc_load)
iproc=proc_dist(iterm)
jproc=proc_dist(jterm)
!Swap terms:
call gaussfit_mpi_add_item(jterm,proc_load(iproc))
call gaussfit_mpi_remove_item(iterm,proc_load(iproc))
call gaussfit_mpi_add_item(iterm,proc_load(jproc))
call gaussfit_mpi_remove_item(jterm,proc_load(jproc))
!Calculate final state
call gaussfit_mpi_calc_deviation(deviation2,nproc,proc_load)
!Swap them only if final state is better than the initial one
if(deviation2<deviation1) then
proc_dist(iterm)=jproc
proc_dist(jterm)=iproc
else
! Return work load to initial state
call gaussfit_mpi_add_item(iterm,proc_load(iproc))
call gaussfit_mpi_remove_item(jterm,proc_load(iproc))
call gaussfit_mpi_add_item(jterm,proc_load(jproc))
call gaussfit_mpi_remove_item(iterm,proc_load(jproc))
! write(*,*)'Back initial state'
! write(*,*)'proc_load',proc_load(:)
end if
end subroutine gaussfit_mpi_swap
!!***
!----------------------------------------------------------------------
!!****f* m_paw_gaussfit/gaussfit_mpi_assign
!! NAME
!! gaussfit_mpi_assign
!!
!! FUNCTION
!! Set task to a processor
!!
!! INPUTS
!!
!! OUTPUT
!!
!! SOURCE
subroutine gaussfit_mpi_assign(iterm,nproc,nterm_bounds,&
& proc_dist,proc_load)
!Arguments ------------------------------------
integer,intent(in)::iterm,nproc,nterm_bounds(2)
integer,intent(inout)::proc_dist(nterm_bounds(1):nterm_bounds(2)),proc_load(nproc)
!Local variables ------------------------------
integer:: iproc,jproc,dev
integer:: deviation(nproc),mindev
character(len=100) :: msg
!************************************************************************
do iproc=1,nproc
!add this term to iproc
call gaussfit_mpi_add_item(iterm,proc_load(iproc))
!calculate the deviation for this configuration
call gaussfit_mpi_calc_deviation(dev,nproc,proc_load)
deviation(iproc)=dev
!remove this term from iproc
call gaussfit_mpi_remove_item(iterm,proc_load(iproc))
end do
!assign to jproc, the proc with minimal deviation above:
jproc=-1; mindev=999999999
do iproc=1,nproc
if(deviation(iproc)<mindev) then
mindev=deviation(iproc)
jproc=iproc
end if
end do
if(jproc==-1) then
! One should not get here!
msg = 'error in accomodate_mpi'
LIBPAW_BUG(msg)
end if
!assign this term for jproc
proc_dist(iterm)=jproc
call gaussfit_mpi_add_item(iterm,proc_load(jproc))
end subroutine gaussfit_mpi_assign
!!***
!----------------------------------------------------------------------
!!****f* m_paw_gaussfit/gaussfit_mpi_main
!! NAME
!! gaussfit_mpi_main
!!
!! FUNCTION
!! Set charge for each processor
!!
!! INPUTS
!!
!! OUTPUT
!!
!! SOURCE
subroutine gaussfit_mpi_main(nproc,nterm_bounds,proc_dist)
!Arguments ------------------------------------
integer,intent(in)::nproc,nterm_bounds(2)
integer,intent(out)::proc_dist(nterm_bounds(1):nterm_bounds(2))
!Local variables ------------------------------
integer:: dev1,dev2,ii
integer:: iproc,iterm,jterm,ngauss,weight
integer:: proc_load(nproc)
character(len=500) :: msg
!************************************************************************
proc_load=0; proc_dist=0 !initializations
!1) get a first-trial distribution:
do iterm=nterm_bounds(2),nterm_bounds(1),-1
call gaussfit_mpi_assign(iterm,nproc,nterm_bounds,proc_dist,proc_load)
end do
!Do the following 20 times
do ii=1,20
! Calculate initial state
call gaussfit_mpi_calc_deviation(dev1,nproc,proc_load)
! Try to swap tasks between two processors:
do iterm=nterm_bounds(2),nterm_bounds(1),-1
do jterm=nterm_bounds(2),nterm_bounds(1),-1
call gaussfit_mpi_swap(iterm,jterm,&
& nproc,nterm_bounds,proc_dist,proc_load)
end do
end do
!! Try to reassign tasks to different processors
do iterm=nterm_bounds(2),nterm_bounds(1),-1
iproc=proc_dist(iterm)
! Remove this job from this node
call gaussfit_mpi_remove_item(iterm,proc_load(iproc))
! Accomodate this again:
call gaussfit_mpi_assign(iterm,nproc,nterm_bounds,proc_dist,proc_load)
end do
! Calculate final state
call gaussfit_mpi_calc_deviation(dev2,nproc,proc_load)
! If final state equals initial state, exit:
if(dev2 == dev1) exit
! write(*,'(a)')'Deviation: ',dev2
end do
!Write down distribution:
write(msg,'(3a)') 'MPI distribution',ch10,'N. gauss, iproc, weight '
call wrtout(std_out,msg,'COLL')
do iterm=nterm_bounds(2),nterm_bounds(1),-1
ngauss=iterm*2
call gaussfit_mpi_set_weight(weight,iterm)
write(msg,'(3(i4,1x))') ngauss,proc_dist(iterm),weight
call wrtout(std_out,msg,'COLL')
end do
write(msg,'(a)') 'Load per processor: '
call wrtout(std_out,msg,'COLL')
do iproc=1,nproc
write(msg,'(i5,1x,i10)') iproc,proc_load(iproc)
call wrtout(std_out,msg,'COLL')
end do
end subroutine gaussfit_mpi_main
!!***
!----------------------------------------------------------------------
!!****f* m_paw_gaussfit/gaussfit_fit
!! NAME
!! gaussfit_fit
!!
!! FUNCTION
!! Fits a given input function f(r) to a sum of gaussians
!!
!! INPUTS
!! chisq= It measures how good is the fitting.
!! It is defined here as sum_i^N_i f(x_i)-y(x_i)/N_i
!! constrains= constraints for Gaussians
!! limit(nparam)= it limits the widths of Gaussians
!! maxiter=maximum number of iterations
!! nparam= number of parameters (a constant times nterm)
!! if(option==1)nparam=nterm*4
!! if(option==2)nparam=nterm*6
!! if(option==3)nparam=nterm*2
!! if(option==4)nparam=nterm*4
!! nterm= number of Gaussians
!! nx=number of points along the x axis
!! option=1 fit to a1 cos(a2 x^2)+ a3 sin( a4 x^2)
!! 2 fit to a1 exp(-a2 x^2)*(a3 cos (a4 x^2) + a5 sin (a6 x^2) )
!! 3 fit to a1 cos (k x^2) + a2 sin (k x^2)
!! 4 fit to a1 exp(-a2 x^2)* (a3 cos(k x^2)+ a4 sin (k x^2))
!! outfile= filename for output (only written if verbosity>1)
!! verbosity= controls output volume
!! weight(nparam)= weights for the fitting procedure
!! x(nx)= points along the x axis
!! y(nx)= function to be fitted
!! rpaw ,optional= paw radius
!!
!! OUTPUT
!! y_out(nx)= fitted function
!!
!! SIDE EFFECTS
!! if(verbosity>1) output files are written with y(x) and y_out(x)
!!
!! SOURCE
subroutine gaussfit_fit(chisq,constrains,&
& limit,maxiter,nparam,nterm,nx,option,outfile,param,&
& verbosity,weight,x,y,y_out)
!Arguments ------------------------------------
integer, intent(in) :: maxiter,nparam
integer,intent(in) :: nterm,nx,option,verbosity
!real(dp),optional,intent(in)::rpaw
!arrays
integer,intent(in)::constrains(nparam)
real(dp),intent(in)::limit(nparam),weight(nparam)
real(dp),intent(in)::x(nx),y(nx)
real(dp),intent(inout)::param(nparam)
real(dp),intent(out)::chisq,y_out(nx)
character(80),intent(in)::outfile
!Local variables ------------------------------
integer, parameter :: wfn_unit=1007
integer::ix
real(dp)::rerror
real(dp),allocatable::sy(:)
! *************************************************************************
LIBPAW_ALLOCATE(sy,(nx))
sy(:)=1.0d0
call gaussfit_rlsf(&
& chisq,constrains,limit,maxiter,&
& nterm,nparam,nx,option,param(1:nparam),&
& verbosity,weight,x,y)
!
if(verbosity>=1) then
if(option==1) then
call gaussfit_calc_deriv_c2(nparam,nterm,nx,1,param,x,y_out)
elseif(option==2) then
call gaussfit_calc_deriv_c(nparam,nterm,nx,1,param,x,y_out)
elseif(option==3) then
call gaussfit_calc_deriv_c3(nparam,nterm,nx,1,param,x,y_out)
elseif(option==4) then
call gaussfit_calc_deriv_c4(nparam,nterm,nx,1,param,x,y_out)
end if
!
!
open(wfn_unit,file=outfile,form='formatted',status='unknown')
! per_error=0.d0
do ix=1, nx
rerror=abs(y(ix)-y_out(ix))
write(wfn_unit,'(6(e20.12,1x))')x(ix),y(ix),y_out(ix),rerror
end do
close(wfn_unit)
end if
LIBPAW_DEALLOCATE(sy)
end subroutine gaussfit_fit
!!***
!----------------------------------------------------------------------
!!****f* m_paw_gaussfit/gaussfit_calc_deriv_r
!! NAME
!! gaussfit_calc_deriv_r
!!
!! FUNCTION
!! Calculate derivatives for Gaussians
!! Only relevant for fitting Gaussians algorithm.
!! The Gaussians expressions are defined in the comments of "gaussfit_main"
!!
!! INPUTS
!!
!! OUTPUT
!!
!! SOURCE
subroutine gaussfit_calc_deriv_r(nterm,nparam,nx,opt,param,x,y_out,&
& deriv) ! optional
!Arguments -------------------------------
integer,intent(in)::nx !number of point in the x grid
integer,intent(in)::nparam !number of parameters
integer,intent(in)::nterm !number of gaussian expressions
integer,intent(in)::opt !option:
!1) calculate only f(x)
!2) calculate f(x) and its derivatives
real(dp),intent(in)::param(nparam) !parameters
real(dp),intent(in)::x(nx) !xgrid
real(dp),intent(out)::y_out(nx) !f(x)
real(dp),optional,intent(out)::deriv(nx,nparam) !derivatives
!Local variables-------------------------------
integer::iexp,ii
real(dp)::alpha1(nterm),alpha2(nterm),alpha3(nterm)
real(dp)::term1(nx,nterm)
real(dp)::aux1(nx)
!real(dp)::step
! *********************************************************************
!
!Initialize
!
y_out(:)=0.d0
!
!Get parameters from parameters array:
!
alpha1(:)=param(1:nterm)
alpha2(:)=param(nterm+1:2*nterm)
alpha3(:)=param(2*nterm+1:3*nterm)
!
!alpha3
!set to constant values of x
!step=rpaw/real(nterm,dp)
!do ii=1,nterm
!raux=step*real(ii,dp)
!alpha3(ii)=raux
!end do
!
!
!
!calculate useful quantities
!
do iexp=1,nterm
aux1(:)=-alpha2(iexp)*(x(:)-alpha3(iexp))**2
term1(:,iexp)=alpha1(iexp)*exp(aux1(:))
end do
!
do iexp=1,nterm
y_out(:)=y_out(:)+term1(:,iexp)
end do
!
!Calculate derivatives:
!
if(opt==2) then
!
! alpha1
!
do iexp=1,nterm
aux1(:)=term1(:,iexp)/alpha1(iexp)
deriv(:,iexp)=aux1(:)
end do
!
! alpha2
!
do iexp=1,nterm
ii=nterm+iexp
aux1(:)=-term1(:,iexp)*(x(:)-alpha3(iexp))
deriv(:,ii)=aux1(:)
deriv(:,ii)=0.1d0
end do
!
! alpha3
!
do iexp=1,nterm
ii=2*nterm+iexp
aux1(:)=term1(:,iexp)*2.d0*alpha2(iexp)
aux1(:)=aux1(:)*(x(:)-alpha3(iexp))
deriv(:,ii)=aux1(:)
end do
end if
end subroutine gaussfit_calc_deriv_r
!!***
!----------------------------------------------------------------------
!!****f* m_paw_gaussfit/gaussfit_calc_deriv_c3
!! NAME
!! gaussfit_calc_deriv_c3
!!
!! FUNCTION
!! Calculate expressions and derivatives for Gaussians fitting.
!! The Gaussians expressions are defined in the comments of "gaussfit_main"
!!
!! INPUTS
!!
!! OUTPUT
!!
!! SOURCE
subroutine gaussfit_calc_deriv_c3(nparam,nterm,nx,opt,param,x,y_out,&
& deriv) ! optional
!Arguments -------------------------------
integer,intent(in)::nparam !number of parameters
integer,intent(in)::nterm !number of gaussian expressions
integer,intent(in)::nx !number of point in the x grid
integer,intent(in)::opt !option:
!1) calculate only f(x)
!2) calculate f(x) and its derivatives
real(dp),intent(in)::param(nparam) !parameters
real(dp),intent(in)::x(nx) !xgrid
real(dp),intent(out)::y_out(nx) !f(x)
real(dp),optional,intent(out)::deriv(nx,nparam) !derivatives
!Local variables-------------------------------
integer::iexp,ii
real(dp)::sep
real(dp)::alpha1(nterm),alpha2(nterm),alpha3(nterm)
real(dp)::term1(nx,nterm)
real(dp)::sin1(nx,nterm),cos1(nx,nterm)
real(dp)::aux1(nx),aux2(nx)
! *********************************************************************
!
!Initialize
!
y_out(:)=0.d0
!
sep=1.2d0
!
!Get param from param array:
!
alpha1(:)=param(1:nterm)
alpha2(:)=param(nterm+1:2*nterm)
!
do ii=1,nterm
alpha3(ii)=sep**(ii)
end do
!
!calculate useful quantities
!
do iexp=1,nterm
aux1(:)=alpha3(iexp)*x(:)**2
!
sin1(:,iexp)=sin(aux1(:))
cos1(:,iexp)=cos(aux1(:))
end do
!
do iexp=1,nterm
aux1(:)=alpha1(iexp)*sin1(:,iexp)
aux2(:)=alpha2(iexp)*cos1(:,iexp)
term1(:,iexp)=aux1(:)+aux2(:)
y_out(:)=y_out(:)+term1(:,iexp)
end do
!
!Calculate derivatives:
!
if(opt==2) then