-
Notifications
You must be signed in to change notification settings - Fork 0
/
spline_cof.f90
1021 lines (900 loc) · 30.6 KB
/
spline_cof.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
!***********************************************************************
!
! routines for calculating spline coefficients
! drivers
!
! Author: Bernhard Seiwald
! Date: 16.12.2000
! 05.11.2001
!
!***********************************************************************
!***********************************************************************
!
! routines for third order spline
!
!***********************************************************************
! ------ third order spline: with testfunction, LSQ, smoothing
!
! AUTHOR: Bernhard Seiwald
!
! DATE: 05.07.2001
SUBROUTINE splinecof3_a(x, y, c1, cn, lambda1, indx, sw1, sw2, &
a, b, c, d, m, f)
!-----------------------------------------------------------------------
!
! compute coefs for smoothing spline with leading function f(x)
! positions of intervals are given by indx
!
! if dabs(c1) > 1e30 -> c1 = 0.0D0
! if dabs(cn) > 1e30 -> cn = 0.0D0
!
! INPUT:
! INTEGER(I4B) , DIMENSION(len_indx) :: indx ... index vector
! contains index of grid points
! ATTENTION:
! x(1),y(1) and x(len_x),y(len_x)
! must be gridpoints!!!
! REAL (kind=dp), DIMENSION(len_x) :: x ...... x values
! REAL (kind=dp), DIMENSION(len_x) :: y ...... y values
! REAL (kind=dp) :: c1, cn .... 1. and last 2. derivative
! REAL (kind=dp), DIMENSION(len_indx) :: lambda . weight for 3. derivative
! INTEGER(I4B) :: sw1 ....
! = 1 -> c1 = 1. deriv 1. point
! = 2 -> c1 = 2. deriv 1. point
! = 3 -> c1 = 1. deriv N. point
! = 4 -> c1 = 2. deriv N. point
! INTEGER(I4B) :: sw2 ....
! = 1 -> cn = 1. deriv 1. point
! = 2 -> cn = 2. deriv 1. point
! = 3 -> cn = 1. deriv N. point
! = 4 -> cn = 2. deriv N. point
! REAL (kind=dp) :: m ...... powers of leading term
! REAL (kind=dp) :: f ...... test function
!
! OUTPUT:
! REAL (kind=dp), DIMENSION(len_indx) :: a, b, c, d ... spline coefs
!
! INTERNAL:
! INTEGER(I4B), PARAMETER :: VAR = 7 ... no of variables
!
! NEEDS:
! solve_systems, calc_opt_lambda3
!
!-----------------------------------------------------------------------
!
!-----------------------------------------------------------------------
! Modules
!-----------------------------------------------------------------------
USE inter_precision, ONLY: I4B, DP
USE inter_interfaces, ONLY: calc_opt_lambda3
USE solve_systems
USE sparse_mod, ONLY: sparse_solve
!-----------------------------------------------------------------------
IMPLICIT NONE
REAL(DP), INTENT(INOUT) :: c1, cn
REAL(DP), DIMENSION(:), INTENT(IN) :: x
REAL(DP), DIMENSION(:), INTENT(IN) :: y
REAL(DP), DIMENSION(:), INTENT(IN) :: lambda1
INTEGER(I4B), DIMENSION(:), INTENT(IN) :: indx
REAL(DP), DIMENSION(:), INTENT(OUT) :: a, b, c, d
INTEGER(I4B), INTENT(IN) :: sw1, sw2
REAL(DP), INTENT(IN) :: m
INTERFACE
FUNCTION f(x,m)
USE inter_precision, ONLY: DP
IMPLICIT NONE
REAL(DP), INTENT(IN) :: x
REAL(DP), INTENT(IN) :: m
REAL(DP) :: f
END FUNCTION f
END INTERFACE
INTEGER(I4B), PARAMETER :: VAR = 7
INTEGER(I4B) :: dim
INTEGER(I4B) :: i_alloc, info
INTEGER(I4B) :: len_x, len_indx
INTEGER(I4B) :: i, j, l, ii, ie
INTEGER(I4B) :: mu1, mu2, nu1, nu2
INTEGER(I4B) :: sig1, sig2, rho1, rho2
INTEGER(I4B), DIMENSION(:), ALLOCATABLE :: indx_lu
REAL(DP) :: h, h_j, x_h, help_i, help_inh
REAL(DP) :: help_a, help_b, help_c, help_d
REAL(DP), DIMENSION(:,:), ALLOCATABLE :: MA
REAL(DP), DIMENSION(:), ALLOCATABLE :: inh, simqa, lambda, omega
len_x = SIZE(x)
len_indx = SIZE(indx)
dim = VAR * len_indx - 2
ALLOCATE(MA(dim, dim), stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3: Allocation for arrays 1 failed!'
ALLOCATE(inh(dim), indx_lu(dim), stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3: Allocation for arrays 2 failed!'
ALLOCATE(simqa(dim*dim), stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3: Allocation for arrays 3 failed!'
ALLOCATE(lambda(SIZE(lambda1)), stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3: Allocation for lambda failed!'
ALLOCATE(omega(SIZE(lambda1)), stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3: Allocation for omega failed!'
!-----------------------------------------------------------------------
IF ( .NOT. ( SIZE(x) == SIZE(y) ) ) THEN
WRITE (*,*) 'splinecof3: assertion 1 failed'
STOP 'program terminated'
END IF
IF ( .NOT. ( SIZE(a) == SIZE(b) .AND. SIZE(a) == SIZE(c) &
.AND. SIZE(a) == SIZE(d) .AND. SIZE(a) == SIZE(indx) &
.AND. SIZE(a) == SIZE(lambda1) ) ) THEN
WRITE (*,*) 'splinecof3: assertion 2 failed'
STOP 'program terminated'
END IF
IF (DABS(c1) > 1.0E30) THEN
c1 = 0.0D0;
END IF
IF (DABS(cn) > 1.0E30) THEN
cn = 0.0D0;
END IF
! setting all to zero
MA(:,:) = 0.0D0
inh(:) = 0.0D0
! check whether points are monotonously increasing or not
DO i = 1, len_x-1
IF (x(i) >= x(i+1)) THEN
PRINT *, 'SPLINECOF3: error i, x(i), x(i+1)', &
i, x(i), x(i+1)
STOP 'SPLINECOF3: error wrong order of x(i)'
END IF
END DO
! check indx
DO i = 1, len_indx-1
IF (indx(i) < 1) THEN
PRINT *, 'SPLINECOF3: error i, indx(i)', i, indx(i)
STOP 'SPLINECOF3: error indx(i) < 1'
END IF
IF (indx(i) >= indx(i+1)) THEN
PRINT *, 'SPLINECOF3: error i, indx(i), indx(i+1)', &
i, indx(i), indx(i+1)
STOP 'SPLINECOF3: error wrong order of indx(i)'
END IF
IF (indx(i) > len_x) THEN
PRINT *, 'SPLINECOF3: error i, indx(i), indx(i+1)', &
i, indx(i), indx(i+1)
STOP 'SPLINECOF3: error indx(i) > len_x'
END IF
END DO
IF (indx(len_indx) < 1) THEN
PRINT *, 'SPLINECOF3: error len_indx, indx(len_indx)', &
len_indx, indx(len_indx)
STOP 'SPLINECOF3: error indx(max) < 1'
END IF
IF (indx(len_indx) > len_x) THEN
PRINT *, 'SPLINECOF3: error len_indx, indx(len_indx)', &
len_indx, indx(len_indx)
STOP 'SPLINECOF3: error indx(max) > len_x'
END IF
! calculate optimal weights for smooting (lambda)
IF ( MAXVAL(lambda1) < 0.0D0 ) THEN
CALL calc_opt_lambda3(x, y, omega)
ELSE
omega = lambda1
END IF
lambda = 1.0D0 - omega
IF (sw1 == sw2) THEN
STOP 'SPLINECOF3: error two identical boundary conditions'
END IF
IF (sw1 == 1) THEN
mu1 = 1
nu1 = 0
sig1 = 0
rho1 = 0
ELSE IF (sw1 == 2) THEN
mu1 = 0
nu1 = 1
sig1 = 0
rho1 = 0
ELSE IF (sw1 == 3) THEN
mu1 = 0
nu1 = 0
sig1 = 1
rho1 = 0
ELSE IF (sw1 == 4) THEN
mu1 = 0
nu1 = 0
sig1 = 0
rho1 = 1
ELSE
STOP 'SPLINECOF3: error in using boundary condition 1'
END IF
IF (sw2 == 1) THEN
mu2 = 1
nu2 = 0
sig2 = 0
rho2 = 0
ELSE IF (sw2 == 2) THEN
mu2 = 0
nu2 = 1
sig2 = 0
rho2 = 0
ELSE IF (sw2 == 3) THEN
mu2 = 0
nu2 = 0
sig2 = 1
rho2 = 0
ELSE IF (sw2 == 4) THEN
mu2 = 0
nu2 = 0
sig2 = 0
rho2 = 1
ELSE
STOP 'SPLINECOF3: error in using boundary condition 2'
END IF
! coefs for first point
i = 0
j = 1
ii = indx((j-1)/VAR+1)
ie = indx((j-1)/VAR+2) - 1
h = x(indx((j-1)/VAR+2)) - x(ii)
! boundary condition 1
i = i + 1
MA(i, 2) = DBLE(mu1)
MA(i, 3) = DBLE(nu1)
MA(i, (len_indx-1)*VAR + 2) = DBLE(sig1)
MA(i, (len_indx-1)*VAR + 3) = DBLE(rho1)
inh(i) = c1
! A_i
i = i + 1
MA(i, j+0 +0) = 1.0D0
MA(i, j+0 +1) = h
MA(i, j+0 +2) = h * h
MA(i, j+0 +3) = h * h * h
MA(i, j+VAR+0) = -1.0D0
! B_i
i = i + 1
MA(i, j+0 +1) = 1.0D0
MA(i, j+0 +2) = 2.0D0 * h
MA(i, j+0 +3) = 3.0D0 * h * h
MA(i, j+VAR+1) = -1.0D0
! C_i
i = i + 1
MA(i, j+0 +2) = 1.0D0
MA(i, j+0 +3) = 3.0D0 * h
MA(i, j+VAR+2) = -1.0D0
! delta a_i
i = i + 1
help_a = 0.0D0
help_b = 0.0D0
help_c = 0.0D0
help_d = 0.0D0
help_i = 0.0D0
DO l = ii, ie
h_j = x(l) - x(ii)
x_h = f(x(l),m) * f(x(l),m)
help_a = help_a + x_h
help_b = help_b + h_j * x_h
help_c = help_c + h_j * h_j * x_h
help_d = help_d + h_j * h_j * h_j * x_h
help_i = help_i + f(x(l),m) * y(l)
END DO ! DO l = ii, ie
MA(i, j+0 +0) = omega((j-1)/VAR+1) * help_a
MA(i, j+0 +1) = omega((j-1)/VAR+1) * help_b
MA(i, j+0 +2) = omega((j-1)/VAR+1) * help_c
MA(i, j+0 +3) = omega((j-1)/VAR+1) * help_d
MA(i, j+0 +4) = 1.0D0
inh(i) = omega((j-1)/VAR+1) * help_i
! delta b_i
i = i + 1
help_a = 0.0D0
help_b = 0.0D0
help_c = 0.0D0
help_d = 0.0D0
help_i = 0.0D0
DO l = ii, ie
h_j = x(l) - x(ii)
x_h = f(x(l),m) * f(x(l),m)
help_a = help_a + h_j * x_h
help_b = help_b + h_j * h_j * x_h
help_c = help_c + h_j * h_j * h_j * x_h
help_d = help_d + h_j * h_j * h_j * h_j * x_h
help_i = help_i + h_j * f(x(l),m) * y(l)
END DO ! DO l = ii, ie
MA(i, j+0 +0) = omega((j-1)/VAR+1) * help_a
MA(i, j+0 +1) = omega((j-1)/VAR+1) * help_b
MA(i, j+0 +2) = omega((j-1)/VAR+1) * help_c
MA(i, j+0 +3) = omega((j-1)/VAR+1) * help_d
MA(i, j+0 +4) = h
MA(i, j+0 +5) = 1.0D0
MA(i, (len_indx-1)*VAR+4) = DBLE(mu1)
MA(i, (len_indx-1)*VAR+5) = DBLE(mu2)
inh(i) = omega((j-1)/VAR+1) * help_i
! delta c_i
i = i + 1
help_a = 0.0D0
help_b = 0.0D0
help_c = 0.0D0
help_d = 0.0D0
help_i = 0.0D0
DO l = ii, ie
h_j = x(l) - x(ii)
x_h = f(x(l),m) * f(x(l),m)
help_a = help_a + h_j * h_j * x_h
help_b = help_b + h_j * h_j * h_j * x_h
help_c = help_c + h_j * h_j * h_j * h_j * x_h
help_d = help_d + h_j * h_j * h_j * h_j * h_j * x_h
help_i = help_i + h_j * h_j * f(x(l),m) * y(l)
END DO ! DO l = ii, ie
MA(i, j+0 +0) = omega((j-1)/VAR+1) * help_a
MA(i, j+0 +1) = omega((j-1)/VAR+1) * help_b
MA(i, j+0 +2) = omega((j-1)/VAR+1) * help_c
MA(i, j+0 +3) = omega((j-1)/VAR+1) * help_d
MA(i, j+0 +4) = h * h
MA(i, j+0 +5) = 2.0D0 * h
MA(i, j+0 +6) = 1.0D0
MA(i, (len_indx-1)*VAR+4) = DBLE(nu1)
MA(i, (len_indx-1)*VAR+5) = DBLE(nu2)
inh(i) = omega((j-1)/VAR+1) * help_i
! delta DELTA d_i
i = i + 1
help_a = 0.0D0
help_b = 0.0D0
help_c = 0.0D0
help_d = 0.0D0
help_i = 0.0D0
DO l = ii, ie
h_j = x(l) - x(ii)
x_h = f(x(l),m) * f(x(l),m)
help_a = help_a + h_j * h_j * h_j * x_h
help_b = help_b + h_j * h_j * h_j * h_j * x_h
help_c = help_c + h_j * h_j * h_j * h_j * h_j * x_h
help_d = help_d + h_j * h_j * h_j * h_j * h_j * h_j * x_h
help_i = help_i + h_j * h_j * h_j * f(x(l),m) * y(l)
END DO ! DO l = ii, ie
MA(i, j+0 +0) = omega((j-1)/VAR+1) * help_a
MA(i, j+0 +1) = omega((j-1)/VAR+1) * help_b
MA(i, j+0 +2) = omega((j-1)/VAR+1) * help_c
MA(i, j+0 +3) = omega((j-1)/VAR+1) * help_d + lambda((j-1)/VAR+1)
MA(i, j+0 +4) = h * h * h
MA(i, j+0 +5) = 3.0D0 * h * h
MA(i, j+0 +6) = 3.0D0 * h
inh(i) = omega((j-1)/VAR+1) * help_i
! coefs for point 2 to len_x_points-1
DO j = VAR+1, VAR*(len_indx-1)-1, VAR
ii = indx((j-1)/VAR+1)
ie = indx((j-1)/VAR+2) - 1
h = x(indx((j-1)/VAR+2)) - x(ii)
! A_i
i = i + 1
MA(i, j+0 +0) = 1.0D0
MA(i, j+0 +1) = h
MA(i, j+0 +2) = h * h
MA(i, j+0 +3) = h * h * h
MA(i, j+VAR+0) = -1.0D0
! B_i
i = i + 1
MA(i, j+0 +1) = 1.0D0
MA(i, j+0 +2) = 2.0D0 * h
MA(i, j+0 +3) = 3.0D0 * h * h
MA(i, j+VAR+1) = -1.0D0
! C_i
i = i + 1
MA(i, j+0 +2) = 1.0D0
MA(i, j+0 +3) = 3.0D0 * h
MA(i, j+VAR+2) = -1.0D0
! delta a_i
i = i + 1
help_a = 0.0D0
help_b = 0.0D0
help_c = 0.0D0
help_d = 0.0D0
help_i = 0.0D0
DO l = ii, ie
h_j = x(l) - x(ii)
x_h = f(x(l),m) * f(x(l),m)
help_a = help_a + x_h
help_b = help_b + h_j * x_h
help_c = help_c + h_j * h_j * x_h
help_d = help_d + h_j * h_j * h_j * x_h
help_i = help_i + f(x(l),m) * y(l)
END DO ! DO l = ii, ie
MA(i, j+0 +0) = omega((j-1)/VAR+1) * help_a
MA(i, j+0 +1) = omega((j-1)/VAR+1) * help_b
MA(i, j+0 +2) = omega((j-1)/VAR+1) * help_c
MA(i, j+0 +3) = omega((j-1)/VAR+1) * help_d
MA(i, j+0 +4) = 1.0D0
MA(i, j-VAR+4) = -1.0D0
inh(i) = omega((j-1)/VAR+1) * help_i
! delta b_i
i = i + 1
help_a = 0.0D0
help_b = 0.0D0
help_c = 0.0D0
help_d = 0.0D0
help_i = 0.0D0
DO l = ii, ie
h_j = x(l) - x(ii)
x_h = f(x(l),m) * f(x(l),m)
help_a = help_a + h_j * x_h
help_b = help_b + h_j * h_j * x_h
help_c = help_c + h_j * h_j * h_j * x_h
help_d = help_d + h_j * h_j * h_j * h_j * x_h
help_i = help_i + h_j * f(x(l),m) * y(l)
END DO ! DO l = ii, ie
MA(i, j+0 +0) = omega((j-1)/VAR+1) * help_a
MA(i, j+0 +1) = omega((j-1)/VAR+1) * help_b
MA(i, j+0 +2) = omega((j-1)/VAR+1) * help_c
MA(i, j+0 +3) = omega((j-1)/VAR+1) * help_d
MA(i, j+0 +4) = h
MA(i, j+0 +5) = 1.0D0
MA(i, j-VAR+5) = -1.0D0
inh(i) = omega((j-1)/VAR+1) * help_i
! delta c_i
i = i + 1
help_a = 0.0D0
help_b = 0.0D0
help_c = 0.0D0
help_d = 0.0D0
help_i = 0.0D0
DO l = ii, ie
h_j = x(l) - x(ii)
x_h = f(x(l),m) * f(x(l),m)
help_a = help_a + h_j * h_j * x_h
help_b = help_b + h_j * h_j * h_j * x_h
help_c = help_c + h_j * h_j * h_j * h_j * x_h
help_d = help_d + h_j * h_j * h_j * h_j * h_j * x_h
help_i = help_i + h_j * h_j * f(x(l),m) * y(l)
END DO ! DO l = ii, ie
MA(i, j+0 +0) = omega((j-1)/VAR+1) * help_a
MA(i, j+0 +1) = omega((j-1)/VAR+1) * help_b
MA(i, j+0 +2) = omega((j-1)/VAR+1) * help_c
MA(i, j+0 +3) = omega((j-1)/VAR+1) * help_d
MA(i, j+0 +4) = h * h
MA(i, j+0 +5) = 2.0D0 * h
MA(i, j+0 +6) = 1.0D0
MA(i, j-VAR+6) = -1.0D0
inh(i) = omega((j-1)/VAR+1) * help_i
! delta DELTA d_i
i = i + 1
help_a = 0.0D0
help_b = 0.0D0
help_c = 0.0D0
help_d = 0.0D0
help_i = 0.0D0
DO l = ii, ie
h_j = x(l) - x(ii)
x_h = f(x(l),m) * f(x(l),m)
help_a = help_a + h_j * h_j * h_j * x_h
help_b = help_b + h_j * h_j * h_j * h_j * x_h
help_c = help_c + h_j * h_j * h_j * h_j * h_j * x_h
help_d = help_d + h_j * h_j * h_j * h_j * h_j * h_j * x_h
help_i = help_i + h_j * h_j * h_j * f(x(l),m) * y(l)
END DO ! DO l = ii, ie
MA(i, j+0 +0) = omega((j-1)/VAR+1) * help_a
MA(i, j+0 +1) = omega((j-1)/VAR+1) * help_b
MA(i, j+0 +2) = omega((j-1)/VAR+1) * help_c
MA(i, j+0 +3) = omega((j-1)/VAR+1) * help_d + lambda((j-1)/VAR+1)
MA(i, j+0 +4) = h * h * h
MA(i, j+0 +5) = 3.0D0 * h * h
MA(i, j+0 +6) = 3.0D0 * h
inh(i) = omega((j-1)/VAR+1) * help_i
END DO ! DO j = VAR+1, VAR*(len_indx-1)-1, VAR
! last point
! delta a_i
i = i + 1
ii = indx((j-1)/VAR+1)
ie = ii
help_a = 0.0D0
help_inh = 0.0D0
l = ii
help_a = help_a + f(x(l),m) * f(x(l),m)
help_inh = help_inh + f(x(l),m) * y(l)
!!$ MA(i, (len_indx-1)*VAR+1) = help_a ! allowed, if omega != 0
!!$ MA(i, (len_indx-2)*VAR+5) = -1.0D0
!!$ inh(i) = help_inh
MA(i, (len_indx-1)*VAR+1) = omega((j-1)/VAR+1) * help_a
MA(i, (len_indx-2)*VAR+5) = omega((j-1)/VAR+1) * (-1.0D0)
inh(i) = omega((j-1)/VAR+1) * help_inh
! delta b_i
i = i + 1
MA(i, (len_indx-2)*VAR+6) = -1.0D0
MA(i, (len_indx-1)*VAR+4) = DBLE(sig1)
MA(i, (len_indx-1)*VAR+5) = DBLE(sig2)
! delta c_i
i = i + 1
MA(i, (len_indx-2)*VAR+7) = -1.0D0
MA(i, (len_indx-1)*VAR+4) = DBLE(rho1)
MA(i, (len_indx-1)*VAR+5) = DBLE(rho2)
!!$! boundary condition 1
!!$ i = i + 1
!!$ MA(i, 2) = DBLE(mu1)
!!$ MA(i, 3) = DBLE(nu1)
!!$ MA(i, (len_indx-1)*VAR + 2) = DBLE(sig1)
!!$ MA(i, (len_indx-1)*VAR + 3) = DBLE(rho1)
!!$ inh(i) = c1
! boundary condition 2
i = i + 1
MA(i, 2) = DBLE(mu2)
MA(i, 3) = DBLE(nu2)
MA(i, (len_indx-1)*VAR + 2) = DBLE(sig2)
MA(i, (len_indx-1)*VAR + 3) = DBLE(rho2)
inh(i) = cn
! make column vector simqa out of matrix MA
!!$ k = 1
!!$ DO i=1, dim
!!$ DO j=1,dim
!!$ simqa(k) = MA(j,i)
!!$ k = k + 1
!!$ END DO
!!$ END DO
! --- for checking matrix ---
!!$! PRINT *, dim
!!$ DO i=1, dim
!!$! PRINT *, 'i = ',i
!!$ DO j=1,dim
!!$ PRINT *, i, j, MA(i,j)
!!$! PRINT *, MA(i,:)
!!$ END DO
!!$ END DO
!!$! PRINT *,'inhom'
!!$! DO i=1, dim
!!$! PRINT *, i, inh(i)
!!$! PRINT *, inh(i)
!!$! END DO
! ---------------------------
! solve system
! CALL solve_eqsys(MA, inh, info)
! IF (info /= 0) STOP 'splinecof3_a: Singular matrix in solve_eqsys()!'
CALL sparse_solve(MA, inh)
! PRINT *, ' A B C D'
! take a(), b(), c(), d()
DO i = 1, len_indx
a(i) = inh((i-1)*VAR+1)
b(i) = inh((i-1)*VAR+2)
c(i) = inh((i-1)*VAR+3)
d(i) = inh((i-1)*VAR+4)
! PRINT *, a(i), b(i), c(i), d(i)
END DO
DEALLOCATE(MA, stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3: Deallocation for arrays 1 failed!'
DEALLOCATE(inh, indx_lu, stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3: Deallocation for arrays 2 failed!'
DEALLOCATE(simqa, stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3: Deallocation for arrays 3 failed!'
DEALLOCATE(lambda, stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3: Deallocation for lambda failed!'
DEALLOCATE(omega, stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3: Deallocation for omega failed!'
END SUBROUTINE splinecof3_a
SUBROUTINE reconstruction3_a(ai, bi, ci, di, h, a, b, c, d)
!-----------------------------------------------------------------------
!
! reconstruct spline coefficients (a, b, c, d) on x(i)
! h := (x - x_i)
!
! INPUT:
! REAL(DP) :: ai, bi, ci, di ... old coefs
! REAL(DP) :: h ................ h := x(i) - x(i-1)
!
! OUTPUT:
! REAL(DP) :: a, b, c, d ....... new coefs
!
!
!-----------------------------------------------------------------------
!
!-----------------------------------------------------------------------
! Modules
!-----------------------------------------------------------------------
USE inter_precision, ONLY: DP
!-----------------------------------------------------------------------
IMPLICIT NONE
REAL(DP), INTENT(IN) :: ai, bi, ci, di
REAL(DP), INTENT(IN) :: h
REAL(DP), INTENT(OUT) :: a, b, c, d
!-----------------------------------------------------------------------
d = di
c = ci + 3.0D0 * h * di
b = bi + h * (2.0D0 * ci + 3.0D0 * h * di)
a = ai + h * (bi + h * (ci + h * di))
END SUBROUTINE reconstruction3_a
SUBROUTINE splinecof3_lo_driv_a(x, y, c1, cn, lambda, w, indx, &
sw1, sw2, a, b, c, d, m, f)
!-----------------------------------------------------------------------
!
! driver routine for splinecof3 ; used for Rmn, Zmn
!
! INPUT:
! INTEGER(I4B) , DIMENSION(len_indx) :: indx ... index vector
! contains index of grid points
! REAL(DP), DIMENSION(no) :: x ...... x values
! REAL(DP), DIMENSION(no) :: y ...... y values
! REAL(DP) :: c1, cn . 1. and last 2. derivative
! REAL(DP), DIMENSION(ns) :: lambda . weight for 3. derivative
! INTEGER(I4B), DIMENSION(ns) :: w ...... weight for point (0,1)
! INTEGER(I4B) :: sw1 .... = 1 -> c1 = 1. deriv 1. point
! = 2 -> c1 = 2. deriv 1. point
! = 3 -> c1 = 1. deriv N. point
! = 4 -> c1 = 2. deriv N. point
! INTEGER(I4B) :: sw2 .... = 1 -> cn = 1. deriv 1. point
! = 2 -> cn = 2. deriv 1. point
! = 3 -> cn = 1. deriv N. point
! = 4 -> cn = 2. deriv N. point
! REAL(DP) :: m ...... powers of leading term
! REAL(DP) :: f ...... test function
!
! OUTPUT:
! REAL(DP), DIMENSION(ns) :: a ...... spline coefs
! REAL(DP), DIMENSION(ns) :: b ...... spline coefs
! REAL(DP), DIMENSION(ns) :: c ...... spline coefs
! REAL(DP), DIMENSION(ns) :: d ...... spline coefs
!
! INTERNAL:
! INTEGER(I4B), PARAMETER :: VAR = 7 ... no of variables
!
!-----------------------------------------------------------------------
!
!-----------------------------------------------------------------------
! Modules
!-----------------------------------------------------------------------
USE inter_precision, ONLY: I4B, DP
USE inter_interfaces, ONLY: splinecof3, reconstruction3
!-----------------------------------------------------------------------
IMPLICIT NONE
INTEGER(I4B), DIMENSION(:), INTENT(IN) :: indx
REAL(DP), INTENT(IN) :: m
REAL(DP), INTENT(INOUT) :: c1, cn
REAL(DP), DIMENSION(:), INTENT(IN) :: x
REAL(DP), DIMENSION(:), INTENT(IN) :: y
REAL(DP), DIMENSION(:), INTENT(IN) :: lambda
INTEGER(I4B), DIMENSION(:), INTENT(IN) :: w
REAL(DP), DIMENSION(:), INTENT(OUT) :: a, b, c, d
INTEGER(I4B), INTENT(IN) :: sw1, sw2
INTERFACE
FUNCTION f(x,m)
USE inter_precision, ONLY: DP
IMPLICIT NONE
REAL(DP), INTENT(IN) :: x
REAL(DP), INTENT(IN) :: m
REAL(DP) :: f
END FUNCTION f
END INTERFACE
INTEGER(I4B) :: dim, no, ns, len_indx
INTEGER(I4B) :: i, j, ie, i_alloc
INTEGER(I4B) :: shift, shifti, shiftv
INTEGER(I4B), DIMENSION(:), ALLOCATABLE :: hi, indx1
REAL(DP) :: h
REAL(DP), DIMENSION(:), ALLOCATABLE :: xn, yn, lambda1
REAL(DP), DIMENSION(:), ALLOCATABLE :: ai, bi, ci, di
no = SIZE(x)
ns = SIZE(a)
len_indx = SIZE(indx)
!-----------------------------------------------------------------------
dim = SUM(w)
IF (dim == 0) THEN
STOP 'error in splinecof3_lo_driv: w == 0'
END IF
ALLOCATE(ai(dim), bi(dim), ci(dim), di(dim), stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3_lo_driv: allocation for arrays 1 failed!'
ALLOCATE(indx1(dim), lambda1(dim), hi(no), stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3_lo_driv: allocation for arrays 2 failed!'
hi = 1
DO i = 1, SIZE(w)
IF ( (w(i) /= 0) .AND. (w(i) /= 1) ) THEN
STOP 'splinecof3_lo_driv: wrong value for w (0/1)'
END IF
IF ( w(i) == 0 ) THEN
IF ( (i+1) <= SIZE(w) ) THEN
ie = indx(i+1)-1
ELSE
ie = SIZE(hi)
END IF
DO j = indx(i), ie
hi(j) = 0
END DO
END IF
END DO
dim = SUM(hi)
ALLOCATE(xn(dim), yn(dim), stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3_lo_driv: allocation for arrays 3 failed!'
! create new vectors for indx and lambda with respect to skipped points
j = 1
shifti = 0
shiftv = 0
DO i = 1, SIZE(indx)
IF ( j <= SIZE(indx1) ) THEN
indx1(j) = indx(i) - shiftv
lambda1(j) = lambda(i-shifti)
END IF
IF ( w(i) /= 0 ) THEN
j = j + 1
ELSE
shifti = shifti + 1
IF ( i+1 <= SIZE(indx) ) THEN
shiftv = shiftv + indx(i+1) - indx(i)
END IF
END IF
END DO
! create new vectors for x and y with respect to skipped points
j = indx1(1)
DO i = 1, SIZE(hi)
IF ( hi(i) /= 0 ) THEN
xn(j) = x(i)
yn(j) = y(i)
j = j+1
END IF
END DO
CALL splinecof3(xn, yn, c1, cn, lambda1, indx1, sw1, sw2, &
ai, bi, ci, di, m, f)
! find first regular point
shift = 1
DO WHILE ( ( shift <= SIZE(w) ) .AND. ( w(shift) == 0 ) )
shift = shift + 1
END DO
! reconstruct spline coefficients from 0 to first calculated coeff.
IF ( ( shift > 1 ) .AND. ( shift < SIZE(w) ) ) THEN
a(shift) = ai(1)
b(shift) = bi(1)
c(shift) = ci(1)
d(shift) = di(1)
DO i = shift-1, 1, -1
h = x(indx(i)) - x(indx(i+1))
CALL reconstruction3(a(i+1), b(i+1), c(i+1), d(i+1), h, &
a(i), b(i), c(i), d(i))
END DO
END IF
! reconstruct all other spline coefficients if needed
j = 0
DO i = shift, ns
IF (w(i) == 1) THEN
j = j + 1
a(i) = ai(j)
b(i) = bi(j)
c(i) = ci(j)
d(i) = di(j)
ELSE
h = x(indx(i)) - x(indx(i-1))
CALL reconstruction3(a(i-1), b(i-1), c(i-1), d(i-1), h, &
a(i), b(i), c(i), d(i))
END IF
END DO
DEALLOCATE(ai, bi, ci, di, stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3_lo_driv: Deallocation for arrays 1 failed!'
DEALLOCATE(indx1, lambda1, hi, stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3_lo_driv: Deallocation for arrays 2 failed!'
DEALLOCATE(xn, yn, stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3_lo_driv: Deallocation for arrays 3 failed!'
END SUBROUTINE splinecof3_lo_driv_a
SUBROUTINE splinecof3_hi_driv_a(x, y, m, a, b, c, d, indx, f)
!-----------------------------------------------------------------------
!
! driver routine for splinecof3_lo_driv
!
! INPUT:
! INTEGER(I4B) , DIMENSION(len_indx) :: indx ... index vector
! contains index of grid points
! INTEGER(I4B), :: choose_rz 1: calc Rmn; 2: Zmn
! REAL(DP), DIMENSION(no) :: x ...... x values
! REAL(DP), DIMENSION(no,no_cur) :: y ...... y values
! REAL(DP), DIMENSION(no_cur) :: m ...... powers of leading term
! REAL(DP) :: f ...... test function
!
! OUTPUT:
! REAL(DP), DIMENSION(ns,no_cur) :: a ...... spline coefs
! REAL(DP), DIMENSION(ns,no_cur) :: b ...... spline coefs
! REAL(DP), DIMENSION(ns,no_cur) :: c ...... spline coefs
! REAL(DP), DIMENSION(ns,no_cur) :: d ...... spline coefs
! INTERNAL:
! REAL(DP), DIMENSION(ns,no_cur) :: lambda3 . weight for 3. derivative
! INTEGER(I4B), DIMENSION(ns,no_cur) :: w ....... weight for point (0,1)
!
!-----------------------------------------------------------------------
!
!-----------------------------------------------------------------------
! Modules
!-----------------------------------------------------------------------
USE inter_precision, ONLY: I4B, DP
USE inter_interfaces, ONLY: splinecof3_lo_driv
!-----------------------------------------------------------------------
IMPLICIT NONE
INTEGER(I4B), DIMENSION(:), INTENT(IN) :: indx
REAL(DP), DIMENSION(:), INTENT(IN) :: m
REAL(DP), DIMENSION(:), INTENT(IN) :: x
REAL(DP), DIMENSION(:,:), INTENT(IN) :: y
REAL(DP), DIMENSION(:,:), INTENT(OUT) :: a, b, c, d
INTERFACE
FUNCTION f(x,m)
USE inter_precision, ONLY: DP
IMPLICIT NONE
REAL(DP), INTENT(IN) :: x
REAL(DP), INTENT(IN) :: m
REAL(DP) :: f
END FUNCTION f
END INTERFACE
REAL(DP), DIMENSION(:,:), ALLOCATABLE :: lambda3
INTEGER(I4B), DIMENSION(:,:), ALLOCATABLE :: w
INTEGER(I4B) :: ns, no_cur
INTEGER(I4B) :: i, sw1, sw2, i_alloc
REAL(DP) :: c1, cn
!-----------------------------------------------------------------------
ns = SIZE(a,1)
no_cur = SIZE(y,2)
ALLOCATE (lambda3(ns,SIZE(y,2)), w(ns,SIZE(y,2)), stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3_hi_driv: Allocation for arrays failed!'
! lambda3 = -1.0D0 !! automatic smoothing
lambda3 = 1.0D0 !! no smoothing
! weights: w(i)=0/1; if(w(i)==0) ... do not use this point
w = 1
sw1 = 2
sw2 = 4
c1 = 0.0D0
cn = 0.0D0
DO i = 1, no_cur
IF ( m(i) /= 0.0D0 ) THEN
w(1,i) = 0 ! system is not defined at y(0)=0
END IF
CALL splinecof3_lo_driv(x, y(:,i), c1, cn, &
lambda3(:,i), w(:,i), indx, sw1, sw2,&
a(:,i), b(:,i), c(:,i), d(:,i), m(i), f)
END DO
DEALLOCATE (lambda3, w, stat = i_alloc)
IF(i_alloc /= 0) STOP 'splinecof3_hi_driv: Deallocation for arrays failed!'
END SUBROUTINE splinecof3_hi_driv_a
SUBROUTINE calc_opt_lambda3_a(x, y, lambda)
! NO FINAL VERSION NOW!!!!!
! calculate optimal weights for smooting (lambda)
!
!-----------------------------------------------------------------------
! Modules
!-----------------------------------------------------------------------
USE inter_precision, ONLY: I4B, DP
USE inter_interfaces, ONLY: dist_lin
!-----------------------------------------------------------------------
IMPLICIT NONE
REAL(DP), DIMENSION(:), INTENT(IN) :: x, y
REAL(DP), DIMENSION(:), INTENT(OUT) :: lambda
INTEGER(I4B) :: i, no
REAL(DP) :: av_a
REAL(DP) :: ymax, xd(3), yd(3)
!-----------------------------------------------------------------------
no = SIZE(x)
av_a = 0.0D0
ymax = MAXVAL(ABS(y))
IF ( ymax == 0.0D0 ) ymax = 1.0D0
DO i = 1, no
IF ( i == 1 ) THEN
xd(1) = x(2)
xd(2) = x(1)
xd(3) = x(3)
yd(1) = y(2)
yd(2) = y(1)
yd(3) = y(3)
CALL dist_lin(xd, yd, ymax, av_a)
ELSE IF ( i == no ) THEN
xd(1) = x(no-2)
xd(2) = x(no)
xd(3) = x(no-1)
yd(1) = y(no-2)
yd(2) = y(no)
yd(3) = y(no-1)
CALL dist_lin(xd, yd, ymax, av_a)
ELSE
CALL dist_lin(x(i-1:i+1), y(i-1:i+1), ymax, av_a)
END IF
!!$ IF (x(i) < 0.2) THEN
!!$ lambda(i) = 1.0D0 - av_a**(1.5D0)
lambda(i) = 1.0D0 - av_a**3
!!$ lambda(i) = 1.0D0 - av_a**(2.5)
!!$ ELSE
!!$ lambda(i) = 1.0D0 - av_a**(5.5)
!!$ END IF
END DO
av_a = SUM(lambda) / DBLE(SIZE(lambda))
lambda = av_a
lambda(1) = 1.0D0
lambda(no) = 1.0D0
END SUBROUTINE calc_opt_lambda3_a