-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathHe_esferico.f90
1085 lines (837 loc) · 26.1 KB
/
He_esferico.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
! Codigo que calculos los autovalores del Helio esfercio.
! El hamiltoniano es:
! H = h(1) + h(2) + lambda* 1/r_>
! donde h(i) = -1/2 d^2 /dr^2 + 1/r
!
! El codigo calcula los autovalores para el estado singlete (espin antiparalelo) por lo tanto
! las funcion de onda espacial es simetrica. Una vez calculado el hamiltoniano en la base producto
! hace una simetrizacion, se reduce el tamaño de la matriz.
!
! El codigo usa B-splines para calcular el hamiltoniano y luega calcula los
! autovalores y autovectores que elija.
!
! El intervalo de integracion es [0, Rmax].
!
! El codigo va variando la variable eta, por lo tanto devuelve los autovalores como
! funcion de eta. Eta representa la intesidad del potencial Vef
!!!!!!!!!!!!!!!!!!!!!! preprocessing variable !!!!!!!
#ifndef Rmax_
#define Rmax_ 10.d0
#endif
#ifndef lum_
#define lum_ 0
#endif
#ifndef c_
#define c_ 0.d0
#endif
#ifndef gamma_
#define gamma_ 0.d0
#endif
#ifndef kord_
#define kord_ 5
#endif
#ifndef l_
#define l_ 30
#endif
#ifndef intg_
#define intg_ 100
#endif
#ifndef nev_
#define nev_ 10
#endif
#ifndef etai_
#define etai_ 0.1d0
#endif
#ifndef etaf_
#define etaf_ 1.d0
#endif
#ifndef num_puntos_eta_
#define num_puntos_eta_ 10
#endif
#ifndef me_
#define me_ 1.d0
#endif
!!!!!!!!!!!!!!!!!!!!!! MODULES !!!!!!!!!!!!!!!!!!!!!!!
module carga
integer :: kord, lum, intg, nev, num_puntos_eta
real(8) :: Rmin, Rmax, eta, etai, etaf
real(8) :: me, delta, c, gamma, ll
integer :: l
character(1) :: tip, l1, l2
end module carga
module matrices
integer :: nk, nb
real(8), allocatable, dimension(:) :: norma
real(8), allocatable, dimension(:,:) :: s, v01, ke, aux
real(8), allocatable, dimension(:,:,:,:) :: Vef
integer, allocatable, dimension(:,:) :: ist, irt, tri
end module matrices
module integracion
integer, allocatable, dimension(:) :: k
real(8), allocatable, dimension(:) :: t, sp, dsp
real(8), allocatable, dimension(:,:) :: x, w, pl
end module integracion
!!!!!!!!!!!!!!!!!!!!!!!!!!!! MAIN !!!!!!!!!!!!!!!!!!!!!!!!!
program Bsplines
use carga
use matrices
use integracion
implicit none
integer :: i, j
real(8) :: time, t1, t2, t3, t4, t5, t6
real(8) :: omp_get_wtime
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
Rmin = 0.d0
Rmax = 10.0 ! valores de inicio y final del intervalo de integracion
tip = 'u' ! tipo de distribucion de knots, u, e, m
lum = 0 ! # de intervalors u en m
c = 0.d0 ! parametod en m dist u en [a,c] y e en (c,b]
gamma = 0.d0 ! parametro de la exponencial
kord = 5 ! orden de los B-splines
l = 50 ! # de intervalos
intg = 100 ! grado de la intregracion de cuadratura gaussiana, >=k
nev = 10 ! # de autovalores que queremos calculara
etai = 0.1d0; etaf = 1.d0
num_puntos_eta = 10
me = 1.d0
if( intg<kord ) intg = kord
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
nk = l+2*kord-1 ! # de knots
! No incluimos primer y ultimo spline debido a psi(0)=psi(R)=0
nb = kord+l-3 ! size base splines
if( nev<0 ) nev = nb
l1 = char(modulo(l, 10) + 48);
l2 = char(modulo(l, 100)/10 + 48);
!###########################################################
!###########################################################
!###########################################################
!open(9,file='electron-autovalores-L='//l2//l1//'.dat')
open(9,file='autovalores-prueba-2.dat')
!###########################################################
!###########################################################
write(9,*) '# El codigo usa B-splines para calcular el hamiltoniano y'
write(9,*) '# calcula los autovalores y autovectores que elija.'
write(9,*) '# El intervalo de integracion es [0, Rmax].'
write(9,*) '# El codigo va variando la variable eta, devuelve los autovalores como'
write(9,*) '# funcion de eta. Eta representa la intesidad del potencial Vef'
write(9,*) '# tip: entrar el tipo de distribucion, u, e, m'
write(9,*) '# tip =', tip
write(9,*) '# lum: entrar el # de intervalos u en m'
write(9,*) '# lum =', lum
write(9,*) '# c: parametro en m dist u en [a,c] y e en (c,b]'
write(9,*) '# c =', c
write(9,*) '# gamma: parametro en exponencial'
write(9,*) '# gamma =', gamma
write(9,*) '# kord: orden de los b-splines'
write(9,*) '# kord =', kord
write(9,*) '# l: # de intervalos'
write(9,*) '# l =', l
write(9,*) '# nb: tamaño base'
write(9,*) '# nb =', nb
write(9,*) '# intg: grado de integracion de la cuadratura >= k'
write(9,*) '# intg =', intg
write(9,*) '# nev: # de autovalores a calcular <=nb'
write(9,*) '# nev =', nev
write(9,*) '# rango de valores de eta='
write(9,*) '# etai =', etai, 'etaf =', etaf
write(9,*) '# masa de electron me =', me
write(9,*) '# autovalores calculados'
!###########################################################
!###########################################################
!###########################################################
close(9)
allocate( Sp(kord), dsp(kord-1))
allocate( x(l,intg), w(l,intg), pl(l,intg))
allocate( t(nk), k(nk))
allocate( norma(nb), s(nb,nb), v01(nb,nb), ke(nb,nb))
allocate( Vef(nb, nb, nb, nb))
t1 = omp_get_wtime();
call KNOTS_PESOS( kord, tip, gamma, Rmin, Rmax, c, l, lum, intg, t, k, x, w, pl)
t2 = omp_get_wtime();
!do i = 1, nk
! write(*,*) i, t(i), k(i)
!end do
!!write(*,*)"calcula las matrices de una particulas"
t3 = omp_get_wtime();
call matrix_elements( )
t4 = omp_get_wtime();
!!write(*,*)"calcula la interaccion"
call interaccion( )
t5 = omp_get_wtime();
do i = 1,nb
do j = i+1,nb
s(j,i) = s(i,j)
v01(j,i) = v01(i,j)
ke(j,i) = ke(i,j)
end do
end do
!!write(*,*)"diagonaliza"
call init_e( )
t6 = omp_get_wtime();
write(*, *) t6 - t1
!write(*,*) " tiempo de KNOTS_PESOS =", t2-t1
!write(*,*) " tiempo de calculo_matrices =", t3-t2
!write(*,*) " tiempo de calculo_interaccion =", t4-t3
!write(*,*) " tiempo de calculo_autovalores =", t5-t4
!write(*,*) " tiempo total =", t5-t1
deallocate(Sp, dsp, x, w, pl)
deallocate(t, k, norma, s, v01, ke, Vef)
call cpu_time(time)
!write(*,*)time/60.d0
end !termina el main, termina el programa
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine matrix_elements( )
use carga
use matrices
use integracion
implicit none
integer :: i, j, m, n, im, in
real(8) :: bm, bn, rr
s = 0.d0 ; v01 = 0.d0 ; ke = 0.d0
do i = kord, kord+l-1
do j = 1, intg
rr = x(k(i),j)
call bsplvb(t, kord, 1, rr, i, Sp)
do m = 1,kord
im = i-kord+m-1
if( im>0.and.im<nb+1 ) then
do n = 1,kord
in = i-kord+n-1
if( in>0.and.in<nb+1 )then
s(im,in) = s(im,in) + sp(m)*sp(n)*w(k(i),j);
v01(im, in) = v01(im, in) + w(k(i),j)*sp(m)*sp(n)/rr;
endif
end do
endif
end do
end do
end do
do i = kord, kord+l-1
do j = 1, intg
rr = x(k(i),j)
do m = i-kord+1, i
if(m>1.and.m<nb+2)then
do n = m,i
if(n<nb+2)then
call bder(rr, t, bm, bn, kord, nk, m, n, i)
ke(m-1,n-1) = ke(m-1,n-1) + 0.5d0*w(k(i),j)*bm*bn/me
endif
end do
endif
end do
end do
end do
end subroutine matrix_elements
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine interaccion( )
use carga
use matrices
use integracion
implicit none
integer :: i1, i2, j1, j2, i, j,it
integer :: m, n, im, in, mp, imp, np
real(8) :: rr1, rr2, w1, w2
real(8), allocatable, dimension(:,:) :: f, g
allocate(f(nb, nb), g(nb, nb))
Vef = 0.d0;
do i1 = kord, kord+l-1
do j1 = 1, intg
rr1 = x(k(i1), j1); w1 = w(k(i1), j1);
f = 0.d0; g = 0.d0;
do i2 = kord, kord+l-1
do j2 = 1, intg
rr2 = x(k(i2), j2); w2 = w(k(i2),j2);
call bsplvb(t, kord, 1, rr2, i2, Sp)
!write(*,*) rr2
do m = 1, kord
im = i2-kord+m-1
if(im>0 .and. im<nb+1)then
do n = 1, kord
in = i2-kord+n-1
if(in>0 .and. in<nb+1)then
if(rr2.LE.rr1)then
f(im, in) = f(im, in) + Sp(m)*Sp(n)*w2/rr1;
else
g(im, in) = g(im, in) + Sp(m)*Sp(n)*w2/rr2;
end if
end if
end do
end if
end do
end do
end do
Sp = 0.d0;
call bsplvb(t, kord, 1, rr1, i1, Sp)
do m = 1, kord
im = i1-kord+m-1
if(im>0 .and. im<nb+1) then
do mp = 1, kord
imp = i1-kord+mp-1;
if(imp>0 .and. imp<nb+1) then
do n = 1, nb
do np = 1, nb
Vef(im, n, imp, np) = Vef(im, n, imp, np) + Sp(m)*Sp(mp)*w1*(f(n, np) + g(n, np))/dsqrt(s(n,n)*s(np,np));
end do
end do
end if
end do
end if
end do
end do
end do
deallocate(f, g)
end subroutine interaccion
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine init_e( )
use matrices
use integracion
use carga
implicit none
integer(4) :: m, n, np, mp, im
real(8) :: nor, nor1, nor2
do n = 1, nb
do np = 1, nb
nor1 = dsqrt(s(n,n)*s(np,np))
do m = 1, nb
do mp = 1, nb
! nor2 = dsqrt(s(m,m)*s(mp,mp));
Vef(n,m,np,mp) = Vef(n,m,np,mp)/(nor1)!*nor2)
end do
end do
end do
end do
!do m=1, nb
! do im=1, nb
! do n=1, nb
! write(*,*) Vef(m,im,n, :)
! end do
! end do
!end do
!stop
do m = 1, nb
do n = m+1, nb
nor = dsqrt(s(m,m)*s(n,n))
s(m,n) = s(m,n)/nor ; s(n,m) = s(m,n)
ke(m,n) = ke(m,n)/nor ; ke(n,m) = ke(m,n)
v01(m,n) = v01(m,n)/nor ; v01(n,m) = v01(m,n)
end do
ke(m,m) = ke(m,m)/s(m,m)
v01(m,m) = v01(m,m)/s(m,m)
end do
do m = 1, nb
norma(m) = dsqrt(s(m,m))
s(m,m) = 1.d0
end do
call sener( )
return
end subroutine init_e
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! sener.f90 basado en ener.f90
! pero es una subroutine de exp_bs_gs.f90 para calcular el GS
! Calcula E_i de una matriz que se le pasa por mudule
subroutine sener( )
use carga
use matrices
use integracion
implicit none
integer(4) :: i, j, m, n, mp, np, dp
integer(4) :: info, ind, indp, NN
real(8) :: raiz
real(8), allocatable, dimension(:,:) :: auvec, val_exp
real(8), allocatable, dimension(:) :: v, e
real(8), allocatable, dimension(:,:) :: mh, ms, mv, hsim
real(8) :: time, t1, t2, t3, t4, t5, t6
real(8) :: omp_get_wtime
raiz = 1.d0/sqrt(2.d0);
dp = nb
NN = nb*(nb+1)/2
allocate( e(nev), v(nev), auvec( NN,nev), mh(dp,dp), val_exp(nev, nev))
allocate( hsim( NN, NN), ms(NN,NN), mv(NN,NN))
!###########################################################
!###########################################################
!###########################################################
!open(31,file='electron-autovalores-L='//l2//l1//'.dat',position='append')
open(11,file='autovalores-prueba-2.dat',position='append')
open(12,file='expectacion_interaccion.dat')
!###########################################################
!###########################################################
!###########################################################
mh = 0.d0
!!! hamiltoniano de una particula
do m = 1, dp
do n = 1, m
mh(m,n) = ke(m,n) - v01(m,n) ;
mh(n,m) = mh(m,n)
end do
end do
delta = (etaf-etai)/dble(num_puntos_eta)
do i = 0, num_puntos_eta-1
!i = 0
eta = etai + dble(i)*delta;
!! write(*,*) i, eta
hsim = 0.d0; ms = 0.d0;
!!! halmiltoniano de dos particulas ya simetrizado
ind = 1
do n = 1, dp
do m = n, dp
indp = 1
do np = 1, dp
do mp = np, dp
if(m.eq.n .and. mp.eq.np)then
hsim( ind, indp) = 2.d0*s(n,np)*mh(n,np) + dble(eta)*Vef(n,n,np,np);
ms( ind, indp) = s(n,np)*s(n,np);
mv( ind, indp) = Vef(n,n,np,np);
elseif(m.ne.n .and. mp.eq.np )then
hsim( ind, indp) = raiz*( 2.d0*s(m,np)*mh(n,np) + 2.d0*s(n,np)*mh(m,np) &
& + dble(eta)*Vef(n,m,np,np) + dble(eta)*Vef(m,n,np,np) );
ms( ind, indp) = 2.d0*raiz*s(n,np)*s(m,np);
mv( ind, indp) = 2.d0*raiz*(Vef(n,m,np,np) + Vef(m,n,np,np));
elseif(m.eq.n .and. mp.ne.np)then
hsim( ind, indp) = raiz*( 2.d0*s(n,mp)*mh(n,np) + 2.d0*s(n,np)*mh(n,mp) &
& + dble(eta)*Vef(n,n,np,mp) + dble(eta)*Vef(n,n,mp,np) );
ms( ind, indp) = 2.d0*raiz*s(n,np)*s(n,mp);
mv( ind, indp) = 2.d0*raiz*(Vef(n,n,np,mp) + Vef(n,n,mp,np));
else
hsim( ind, indp) = s(n,np)*mh(m,mp) + s(n,mp)*mh(m,np) &
&+ s(m,mp)*mh(n,np) + s(m,np)*mh(n,mp) &
&+ dble(eta)*0.5d0*(Vef(n,m,np,mp) + Vef(n,m,mp,np) + Vef(m,n,np,mp) + Vef(m,n,mp,np));
ms( ind, indp) = s(n,np)*s(m,mp) + s(n,mp)*s(m,np);
mv( ind, indp) = 0.5d0*(Vef(n,m,np,mp) + Vef(n,m,mp,np) + Vef(m,n,np,mp) + Vef(m,n,mp,np));
endif
indp = indp + 1
end do
end do
ind = ind + 1
end do
end do
if( NN.ne.(ind-1) ) stop
if( NN.ne.(indp-1) ) stop
!write(*,*) NN, NN
!do n = 1, NN
! write(*,*) hsim(n, :)!-hsim(m,n)
!end do
!stop;
!t1 = omp_get_wtime();
!call KNOTS_PESOS( kord, tip, gamma, Rmin, Rmax, c, l, lum, intg, t, k, x, w, pl)
call eigen_value( NN, nev, info, hsim, ms, e, auvec)
!t2 = omp_get_wtime();
val_exp = MATMUL( TRANSPOSE(auvec), MATMUL(mv, auvec))
do j = 1, nev
v(j) = val_exp(j,j);
end do
!muestro los valores
write(11,6) eta, (e(m), m = 1, 25)
write(*,*) v
write(*,6) eta, (e(m), m = 1, 2)
!e = 0.d0
!info = 0
end do
close(11)
close(12)
deallocate(e, v, auvec, mh, ms, hsim, mv, val_exp)
6 format(e22.14,1x,1000(1x,e22.14))
return
end subroutine sener
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine eigen_value(dp,nev,info,nh,s,x,avec)
!
! calcula usando LAPACK N=nev autovalores
implicit none
integer*4 dp,nev
real(8),dimension (nev)::x
real(8),dimension (dp,nev)::avec
real(8),dimension (dp,dp)::nh,s
! definiciones para DSYGVX
integer(4)::itype,numau,lwork,info,il,iu
character(1)::jobz,range,uplo
real(8)::vl,vu,abstol
real(8),dimension(dp,dp)::z
real(8),dimension(9*dp)::work
integer(4),dimension(5*dp)::iwork
integer(4),dimension(dp)::ifail
lwork=9*dp ! lwork=dim(work)
!!!!!!!!!!!!!!!!!!!!!!!!!
! reemplaza tred2 y tqli del NR por dsygvx de LAPACK
! parametros para dsygvx:
itype=1 ! especifica que A x=lambda Bx
jobz='V' ! V (N) (NO) calcula autovectores
!range='V'! autovalores en el intervalo [VL,VU]
!vl=-1.01d0 ! pongo vl menor que el menor autovalor
!vu=0.d0 ! por ahora solo interezan los autovalores negativos
uplo='U' ! la matriz es storeada superior
! N es dp
!A es nh
! lda = dp
! B = s ! NO hay que hacer cholesky !!!
! ldb = dp
! il, iu autov entre il y iu SI range='I'
range='I'
il=1
iu=nev
abstol=1.d-12 ! por decir algo (?)aconsejan abstol=2*DLAMCH('S') ; ????
! M = numau# de autovalores calculados : output
! W = x , de menor a mayor
! se debe definir un array Z(LDZ,M) NO conocemos M a priori! poner N=dp
! Z devuelve eigenvectors
! ldz=dp : dimension de Z
! definir work(lwork)
lwork=9*dp
! definir integer array iwork(5*N)
! ifail=output:
!If JOBZ = 'V', then if INFO = 0, the first M elements of
! IFAIL are zero. If INFO > 0, then IFAIL contains the
! indices of the eigenvectors that failed to converge.
! If JOBZ = 'N', then IFAIL is not referenced.
! INFO (output) INTEGER
! = 0: successful exit
! < 0: if INFO = -i, the i-th argument had an illegal value
! > 0: DPOTRF or DSYEVX returned an error code:
! <= N: if INFO = i, DSYEVX failed to converge;
! i eigenvectors failed to converge. Their indices
! are stored in array IFAIL.
! > N: if INFO = N + i, for 1 <= i <= N, then the leading
! minor of order i of B is not positive definite.
! The factorization of B could not be completed and
! no eigenvalues or eigenvectors were computed.
call DSYGVX( ITYPE, JOBZ, RANGE, UPLO, dp, nh, dp, s, dp, VL, VU, IL, IU, ABSTOL, numau, x, Z, dp, WORK, LWORK, IWORK, IFAIL, INFO )
!!!!!!!!!!!!!!!!!!
avec(1:dp,1:nev)=z(1:dp,1:nev)
return
end subroutine eigen_value
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! bajada de http://www.atom.physto.se/~lindroth/comp08/bget.f
subroutine bder(rr,t,dm,dn,kord,np,indexm,indexn,left)
! returns the value of (d/dx) Bspline(kord,index) in rr
! The first Bspline is called spline #1 (i.e. index=1)
! the first knot point is in t(1)
! np= number of knot points (distinct or multiple) including
! the ghost points: N phyical points np=N +2*(kord-1)
implicit none
integer i,left,kord,indexm,indexn,np
real*8 dm,dn,t(np),Sp(kord),rr
dm=0.d0 ; dn=0.d0
! if rr=t(np) then the routine assumes that
! there is kord knotpoints in the last physical point and
! returns Bder.ne.zero if index is np-kord, or np-kord-1
if(rr.gt.t(np)) then
return
endif
if(rr.lt.t(1)) then
return
endif
! do it=1,np
! if(rr.ge.t(it)) left=it
! end do
if(abs(rr-t(np)).lt.1.d-10) then
!if(index.lt.np-kord-1) return
if(indexm.eq.np-kord) then
dm=dble(kord-1)/(t(np)-t(np-kord))
else if(indexm.eq.np-kord-1) then
dm=-dble(kord-1)/(t(np)-t(np-kord))
end if
if(indexn.eq.np-kord) then
dn=dble(kord-1)/(t(np)-t(np-kord))
else if(indexn.eq.np-kord-1) then
dn=-dble(kord-1)/(t(np)-t(np-kord))
end if
return
end if
call bsplvb(t,kord-1,1,rr,left,Sp)
if(indexm-left+kord.ge.1.or.indexm-left+kord.le.kord)then
i=indexm-left+kord
if(i.eq.1) then
dm=dble(kord-1)*(-Sp(i)/(t(indexm+kord)-t(indexm+1)))
else if(i.eq.kord) then
dm=dble(kord-1)*(Sp(i-1)/(t(indexm+kord-1)-t(indexm)))
else
dm=dble(kord-1)*(Sp(i-1)/(t(indexm+kord-1)-t(indexm))- &
& Sp(i)/(t(indexm+kord)-t(indexm+1)))
end if
endif
if(indexn-left+kord.ge.1.or.indexn-left+kord.le.kord)then
i=indexn-left+kord
if(i.eq.1) then
dn=dble(kord-1)*(-Sp(i)/(t(indexn+kord)-t(indexn+1)))
else if(i.eq.kord) then
dn=dble(kord-1)*(Sp(i-1)/(t(indexn+kord-1)-t(indexn)))
else
dn=dble(kord-1)*(Sp(i-1)/(t(indexn+kord-1)-t(indexn))- &
& Sp(i)/(t(indexn+kord)-t(indexn+1)))
end if
endif
return
end
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
SUBROUTINE bsplvb(t,jhigh,index,x,left,biatx)
! INCLUDE 'standa.typ'
! INCLUDE 'spldim.def'
! include 'ndim_inc'
implicit none
integer(4),PARAMETER::JMAX=100
integer index,jhigh,left,i,j,jp1
real*8 t,x,biatx,deltal,deltar,saved,term
DIMENSION biatx(jhigh),t(left+jhigh),deltal(jmax),deltar(jmax)
SAVE deltal,deltar
DATA j/1/
! write(6,*) ' jmax=',jmax
GO TO (10,20),index
10 j = 1
biatx(1) = 1.d0
IF (j .GE. jhigh) GO TO 99
20 CONTINUE
jp1 = j + 1
deltar(j) = t(left+j) - x
deltal(j) = x - t(left+1-j)
! write(6,'(1pd12.4,2(i5,1pd14.6))')
! : x,left+j,t(left+j),left+1-j,t(left+1-j)
! write(6,'(i3,1p3d12.4)') j,deltal(j),deltar(j),
! : abs(deltal(j)-deltar(j))
saved = 0.d0
DO i = 1,j
! write(6,'(2i3,1p3d12.4)') i,j,deltal(jp1-1),deltar(i),
! : abs(deltal(jp1-1)-deltar(i))
term = biatx(i)/(deltar(i) + deltal(jp1-i))
biatx(i) = saved + deltar(i)*term
saved = deltal(jp1-i)*term
END DO
biatx(jp1) = saved
j = jp1
IF (j .LT. jhigh) GO TO 20
99 RETURN
END
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
SUBROUTINE KNOTS_PESOS(kord,tip,gamma,a,b,c,l,lum,intg,t,k,x,w,pl)
!lo único que hace esta subrutina principal es dividir los casos segun se quiera knots uniformes, exponencial, mixto
!cada una de las subrutinas que llama ademas calcula los x y los w, abscisas y pesos para la cuadratura gaussiana, y los polinomios
! INPUT:
! kord: orden de los b-splines
! tip; character(1): 'u' ; 'e' ; 'm': dist uniforme, exp o mixta de knots
! gamma : param dist. e (no usado si tip='u')
! a
! b ; a<b todo es calculado en el intervalo [a,b]
! c solo usado si tip='m' => a<c<b ; dist u en [a,c]; e en (c,b]; c es a u-knot
! NOTA I : dist optima: m; con c~r_0 donde v(r) es apreciablemente no nulo
! l = # total de sub_intervalos en [a,b]
! lum: solo se usa si tip=m : # subint con dist u => l-lem con dist e
! intg = # de puntos que usamos en la integral x cuadratura en c/subintervalo
! OUTPUT:
! t(l+2*kord-1) ! los (l+2*kord-1) knots (contando multiplicidades en a y b)
! k(l+2*kord-1) ! da el # de intervalo j para el i-esimo nodo 1<=k<=l
! x(l,intg),w(l,intg) : posiciones y pesos de los intg puntos en c/intervalo
! pl(l,intg),w(l,intg) : polinomios de Legendre en cada punto
! NOTA II : los p_l los calcula gratis la subr gauleg de int. x cuad. del NR
! y a veces hacen falta=> version modif. de gauleg que los da como output
implicit none
integer(4)::kord,l,lum,intg,i1,j1
character(1)::tip
real(8)::gamma,a,b,c
integer(4),dimension(l+2*kord-1)::k
real(8),dimension(l+2*kord-1)::t
real(8),dimension(l,intg)::x,w,pl
real(8)::rr2
!!!!!!
if(tip=='u')then
call dist_unif(kord,a,b,l,intg,t,k,x,w,pl)
elseif(tip=='e')then
call dist_exp(kord,gamma,a,b,l,intg,t,k,x,w,pl)
elseif(tip=='m')then
call dist_mix(kord,gamma,a,b,c,l,lum,intg,t,k,x,w,pl)
else
!write(*,*)'error 1 en KNOTS_PESOS :',tip,' no corresponde a una distribucion'
stop
endif
return
end SUBROUTINE KNOTS_PESOS
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
SUBROUTINE DIST_UNIF(kord,a,b,l,intg,t,k,x,w,pl)
implicit none
integer(4)::kord,l,intg
real(8)::a,b
integer(4),dimension(l+2*kord-1)::k
real(8),dimension(l+2*kord-1)::t
real(8),dimension(l,intg)::x,w,pl
!!!!!!
integer(4)::i,j,nk
real(8)::ri,rf,dr
real(8),dimension(intg)::vx,vw,vpl
nk=l+2*kord-1 ! # de knots
dr=(b-a)/dfloat(l)
! calcula los puntos y pesos para la cuadratura
x=0.d0;w=0.d0
do i=1,l
ri=a+dfloat(i-1)*dr
rf=ri+dr
call gauleg_pl(ri,rf,vx,vw,vpl,intg)
do j=1,intg
x(i,j)=vx(j)
w(i,j)=vw(j)
pl(i,j)=vpl(j)
end do
end do
t(1)=a
k(1)=1
do i=2,kord
t(i)=t(1)
k(i)=1
end do
do i=kord+1,kord+l
t(i)=t(i-1)+dr
k(i)=k(i-1)+1
end do
do i=kord+l+1,nk
t(i)=t(i-1)
k(i)=k(i-1)
end do
return
end SUBROUTINE DIST_UNIF
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
SUBROUTINE DIST_EXP(kord,gamma,a,b,l,intg,t,k,x,w,pl)
implicit none
integer(4)::kord,l,intg
real(8)::gamma,a,b
integer(4),dimension(l+2*kord-1)::k
real(8),dimension(l+2*kord-1)::t
real(8),dimension(l,intg)::x,w,pl
!!!!!!
integer(4)::i,j,nk
real(8)::ri,rf,dr,ye
real(8),dimension(intg)::vx,vw,vpl
nk=l+2*kord-1 ! # de knots
dr=(b-a)/(dexp(gamma)-1.d0)
ye=gamma/dfloat(l)
! calcula los puntos y pesos para la cuadratura
x=0.d0;w=0.d0
do i=1,l
ri=a+dr*(dexp(ye*dfloat(i-1))-1.d0)
rf=a+dr*(dexp(ye*dfloat(i))-1.d0)
call gauleg_pl(ri,rf,vx,vw,vpl,intg)
do j=1,intg
x(i,j)=vx(j)
w(i,j)=vw(j)
pl(i,j)=vpl(j)
end do
end do
t(1)=a
k(1)=1
do i=2,kord
t(i)=t(1)
k(i)=1
end do
do i=kord+1,kord+l
t(i)=a+dr*(dexp(ye*dfloat(k(i-1)))-1.d0)
k(i)=k(i-1)+1
end do
do i=kord+l+1,nk
t(i)=t(i-1)
k(i)=k(i-1)
end do
return
end SUBROUTINE DIST_EXP
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
SUBROUTINE DIST_MIX(kord,gamma,a,b,c,l,lum,intg,t,k,x,w,pl)
implicit none
integer(4)::kord,l,lum,intg
real(8)::gamma,a,b,c
integer(4),dimension(l+2*kord-1)::k
real(8),dimension(l+2*kord-1)::t
real(8),dimension(l,intg)::x,w,pl
!!!!!!
integer(4)::i,j,le,nk
real(8),dimension(lum,intg)::xu,wu,plu
real(8),dimension(l-lum,intg)::xe,we,ple
integer(4),dimension(lum+2*kord-1)::ku
real(8),dimension(lum+2*kord-1)::tu
integer(4),dimension(l-lum+2*kord-1)::ke
real(8),dimension(l-lum+2*kord-1)::te