-
Notifications
You must be signed in to change notification settings - Fork 6
/
mo_fit.f90
1772 lines (1437 loc) · 55.3 KB
/
mo_fit.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
MODULE mo_fit
! This module contains the numerical recipes routines to fit arbitrary functions to discrete data
! Usage:
! USE mo_fit, ONLY: fitfun, & ! Fit function: wrapper of svdfit & svdvar
! fpoly, fpoly_sp, fpoly_dp, & ! Polynomial base functions
! polyfit, & ! Polynomial fitting routine using Lapack
! polyval ! Calculates polynomial values for an vectorial input
! svdfit, svdvar ! fitting routines + errors
! Literature
! SVDFit
! WH Press, SA Teukolsky, WT Vetterling, BP Flannery,
! Numerical Recipes in Fortran 90 - The Art of Parallel Scientific Computing, 2nd Edition
! Volume 2 of Fortran Numerical Recipes, Cambridge University Press, UK, 1996
! Polyfit
! http://rosettacode.org/wiki/Polynomial_regression#Fortran
! Attention: polyfit uses LAPACK
! Written March 2011, Matthias Cuntz - copied and adapted numerical recipes subroutines
! - fitfun for easy use with mask
! - linfit Model II: geometric mean regression
! Modified Matthias Cuntz, Nov 2011 - sp, dp
! - documentation
! Gregor Schuldt, Sep 2014 - polyval
! Matthias Cuntz, Jul 2015 - pgiFortran
! Matthias Cuntz, Mai 2016 - rm linfit -> use mo_linfit
! Matthias Cuntz, Feb 2020 - make pgiFortran code for all compilers
! License
! -------
! This file is part of the JAMS Fortran package, distributed under the MIT License.
!
! Copyright (c) 2011-2020 Matthias Cuntz - mc (at) macu (dot) de
!
! Permission is hereby granted, free of charge, to any person obtaining a copy
! of this software and associated documentation files (the "Software"), to deal
! in the Software without restriction, including without limitation the rights
! to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
! copies of the Software, and to permit persons to whom the Software is
! furnished to do so, subject to the following conditions:
!
! The above copyright notice and this permission notice shall be included in all
! copies or substantial portions of the Software.
!
! THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
! IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
! FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
! AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
! LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
! OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
! SOFTWARE.
! Note on Numerical Recipes License
! ---------------------------------
! Be aware that some code is under the Numerical Recipes License 3rd
! edition <http://www.nr.com/aboutNR3license.html>
! The Numerical Recipes Personal Single-User License lets you personally
! use Numerical Recipes code ("the code") on any number of computers,
! but only one computer at a time. You are not permitted to allow anyone
! else to access or use the code. You may, under this license, transfer
! precompiled, executable applications incorporating the code to other,
! unlicensed, persons, providing that (i) the application is
! noncommercial (i.e., does not involve the selling or licensing of the
! application for a fee), and (ii) the application was first developed,
! compiled, and successfully run by you, and (iii) the code is bound
! into the application in such a manner that it cannot be accessed as
! individual routines and cannot practicably be unbound and used in
! other programs. That is, under this license, your application user
! must not be able to use Numerical Recipes code as part of a program
! library or "mix and match" workbench.
! Businesses and organizations that purchase the disk or code download,
! and that thus acquire one or more Numerical Recipes Personal
! Single-User Licenses, may permanently assign those licenses, in the
! number acquired, to individual employees. Such an assignment must be
! made before the code is first used and, once made, it is irrevocable
! and can not be transferred.
! If you do not hold a Numerical Recipes License, this code is only for
! informational and educational purposes but cannot be used.
USE mo_kind, ONLY: i4, sp, dp
implicit NONE
PUBLIC :: fitfun ! Wrapper of svdfit and svdvar
PUBLIC :: fpoly ! Routine to fit polynomial with fitfun or svdfit
PUBLIC :: polyfit ! Fit a polynomial (with Lapack)
PUBLIC :: polyval ! Calculates polynomial-values for an vectorial input
PUBLIC :: svdfit ! Parameter fitting with singular value decomposition
PUBLIC :: svdvar ! Variance of fitted parameters with svdfit
public :: fpoly_sp, fpoly_dp ! individual routines for fpoly to pass to fitting functions
! Public interfaces
! ------------------------------------------------------------------
! NAME
! fitfun
! PURPOSE
! Fits a given function func with the parameters a to the input data by using singular value decomposition.
! The function func must give a vector of the base functions for a given x.
! For example func for a polynomial should return the vector [1,x,x**2,x**3,x**4,...].
! From the numerical recipes documentation for svdfit:
! Given a set of N data points x, y with individual standard deviations sig, all arrays of length
! N, use chi^2 minimization to determine the M coefficients a of a function that depends linearly
! on a, y = sum_i=1^M(a x func_i(x)), with x meaning the dot-product.
! Here we solve the fitting equations using singular value decomposition of the NxM matrix, as in paragraph 2.6.
! The program returns the M basis functions evaluated at x=X in fit, the M fit parameters a,
! their sigma-errors siga, and chi^2 chisq.
! If an optinal mask is given, the function is only fit at those locations that correspond to true values in the mask.
! CALLING SEQUENCE
! call fitfun(xin, yin, sigin, a, func, mask=mask, fit=fit, chisq=chisq, siga=siga)
! INTENT(IN)
! real(sp/dp) :: xin(:) 1D-array with input x
! real(sp/dp) :: yin(:) 1D-array with input y
! real(sp/dp) :: sigin(:) 1D-array with input sigma on y
! real(sp/dp), dimension(M) :: FUNCTION func(x,n) Function that outputs the n basis functions
! INTENT(INOUT)
! None
! INTENT(OUT)
! real(sp/dp), dimension(M) :: a fitted M parameters
! INTENT(IN), OPTIONAL
! logical :: mask(:) 1D-array of logical values with size(vec).
! If present, only those locations in vec corresponding to the true values in mask are used.
! INTENT(INOUT), OPTIONAL
! None
! INTENT(OUT), OPTIONAL
! real(sp/dp), dimension(N) :: fit Function evaluated at x=X
! real(sp/dp), dimension(M) :: siga Error of fitted parameters
! real(sp/dp) :: chisq Minimum chi^2
! RESTRICTIONS
! None
! EXAMPLE
! ! fit polynomial with degree 5, then detrend original
! integer(i4), parameter :: npt=100, npol=4 ! npol=degrees+1
! real(dp), dimension(npt) :: x, y, sig
! real(dp), dimension(npol) :: a
! logical, dimension(npt) :: maske
! maske = .true.
! maske(10) = .false.
! ! fit
! call fitfun(x, y, sig, a, fpoly_dp, fit=fit, mask=maske)
! ! detrend
! do i=1, npt
! y(i) = y(i) - dot_product(a, fpoly(x(i),npol))
! end do
! ! or
! y = y - fit
! LITERATURE
! Press WH, Teukolsky SA, Vetterling WT, & Flannery BP - Numerical Recipes in Fortran 90 -
! The Art of Parallel Scientific Computing, 2nd Edition, Volume 2 of Fortran Numerical Recipes,
! Cambridge University Press, UK, 1996
! HISTORY
! Written, Matthias Cuntz, Mar 2011
INTERFACE fitfun
MODULE PROCEDURE fitfun_sp, fitfun_dp
END INTERFACE fitfun
! ------------------------------------------------------------------
! NAME
! fpoly
! PURPOSE
! Fitting routine for a polynomial of degree n-1, i.e. with n coefficients.
! Returns vector with [1,x,x**2,x**3,...,x**(n-1)].
! To use with fitfun or svdfit.
! CALLING SEQUENCE
! vec = func(x,n)
! INTENT(IN)
! real(sp/dp) :: x x
! real(sp/dp) :: n n-1 powers of x
! INTENT(INOUT)
! None
! INTENT(OUT)
! real(sp/dp), dimension(n) :: func vec with [1,x,x**2,x**3,...,x**(n-1)]
! INTENT(IN), OPTIONAL
! None
! INTENT(INOUT), OPTIONAL
! None
! INTENT(OUT), OPTIONAL
! None
! RESTRICTIONS
! None
! EXAMPLE
! -> see example at fitfun
! LITERATURE
! Press WH, Teukolsky SA, Vetterling WT, & Flannery BP - Numerical Recipes in Fortran 90 -
! The Art of Parallel Scientific Computing, 2nd Edition, Volume 2 of Fortran Numerical Recipes,
! Cambridge University Press, UK, 1996
! HISTORY
! Written, Matthias Cuntz, Mar 2011
INTERFACE fpoly
MODULE PROCEDURE fpoly_sp, fpoly_dp
END INTERFACE fpoly
! ------------------------------------------------------------------
! NAME
! polyfit
! PURPOSE
! Find an approximating polynom of known degree for a given data.
! The function polyfit implements least squares approximation of a function
! defined in the points as specified by the arrays xi and yi.
! The basis dj is xj, j=0,1,..,N. The implementation is straightforward.
! First the plane matrix A is created. Aji=dj(xi). Then the linear problem
! AATc=Ay is solved. The result cj are the coefficients. Constraint_Error
! is propagated when dimensions of X and Y differ or else when the problem is ill-defined.
! If an optinal mask is given, the function is only fit at those locations that correspond to true values in the mask.
! CALLING SEQUENCE
! out = polyfit(x, y, d, mask=mask)
! INTENT(IN)
! real(sp/dp) :: x(:) x
! real(sp/dp) :: y(:) y
! integer :: d d+1 parameters = polynom of degree d
! INTENT(INOUT)
! None
! INTENT(OUT)
! real(sp/dp), dimension(d+1) :: polyfit fitted d+1 parameters for x^0, x^1, ..., x^d
! INTENT(IN), OPTIONAL
! logical :: mask(:) 1D-array of logical values with size(x)
! If present, only those locations in vec corresponding to the true values in mask are used.
! INTENT(INOUT), OPTIONAL
! None
! INTENT(OUT), OPTIONAL
! None
! RESTRICTIONS
! None
! EXAMPLE
! integer, parameter :: degree = 2
! integer :: i
! real(dp), dimension(11) :: x = (/ (i,i=0,10) /)
! real(dp), dimension(11) :: y = (/ 1, 6, 17, 34, 57, 86, 121, 162, 209, 262, 321 /)
! real(dp), dimension(degree+1) :: a
! a = polyfit(x, y, degree)
! write (*, '(F9.4)') a
! -> gives 1.0, 2.0, 3.0
! LITERATURE
! http://rosettacode.org/wiki/Polynomial_regression#Fortran
! HISTORY
! Written, Matthias Cuntz, Mar 2011
INTERFACE polyfit
MODULE PROCEDURE polyfit_sp, polyfit_dp
END INTERFACE polyfit
! ------------------------------------------------------------------
! NAME
! polyval
! PURPOSE
! Evaluates a polynomial for a given vector.
! Input 'poly' must be an array containing the values for the coefficients of the polynomial.
! For example 'poly' for the polynomial "1 + 2 * x + 3 * x**2" should be (/ 1, 2, 3 /)
! Input 'xwert' can be scalar or 1D-array with values for which the user wants to evaluate the given polynomial.
! CALLING SEQUENCE
! vec = polyval(poly,xwert)
! INTENT(IN)
! real(sp/dp), dimension(:) :: poly polynomial-coefficients
! real(sp/dp)[, dimension(:)] :: xwert Evaluation-value(s)
! INTENT(INOUT)
! None
! INTENT(OUT)
! real(sp/dp)[, dimension(size(xwert,1))] :: polyval (vec with) polynomial-value(s) at xwert
! INTENT(IN), OPTIONAL
! None
! INTENT(INOUT), OPTIONAL
! None
! INTENT(OUT), OPTIONAL
! None
! RESTRICTIONS
! None
! EXAMPLE
! For the polynom "1 + 2 * x + 3 * x**2" and the values 1 and 2
! poly = (/ 1, 2, 3/)
! xwert = (/ 1, 2 /)
! Erg = polyval(poly,xwert)
! ----> Erg = (/ 6, 17 /)
! LITERATURE
!
! HISTORY
! Written, Gregor Schuldt, Sep 2014
INTERFACE polyval
MODULE PROCEDURE polyval_0d_sp, polyval_1d_sp, polyval_0d_dp, polyval_1d_dp
END INTERFACE polyval
! ------------------------------------------------------------------
! NAME
! svdfit
! PURPOSE
! Fits a given function func with the parameters a to the input data by using singular value decomposition.
! The function func must give a vector of the base functions for a given x.
! For example func for a polynomial should return the vector [1,x,x**2,x**3,x**4,...].
! From the numerical recipes documentation:
! Given a set of datapoints x(1:ndata), y(1:ndata) with individual standard deviations
! sig(1:ndata), use chi2 minimization to determine the ma coefficients a of the fitting function
! y=sum(ai*afunci(x))). Here we solve the fitting equations using singular value decomposition
! of the n data by ma matrix, as in paragraph 2.6. Arrays u(1:mp,1:np), v(1:np,1:np),
! w(1:np) provide workspace on input; on output they define the singular value decomposition,
! and can be used to obtain the covariance matrix. mp,np are the physical dimensions
! of the matrices u,v,w, as indicated above. It is necessary that mp>=ndata, np>=ma. The
! program returns values for the ma fit parameters a, and chi2, chisq. The user supplies a
! subroutine funcs(x,afunc,ma) that returns them a basis functions evaluated at x=x
! in the array afunc.
! CALLING SEQUENCE
! call svdfit(x, y, sig, a, v, w, chisq, func)
! INTENT(IN)
! real(sp/dp) :: x(:) 1D-array with input x
! real(sp/dp) :: y(:) 1D-array with input y
! real(sp/dp) :: sig(:) 1D-array with input sigma on y
! real(sp/dp), dimension(M) :: FUNCTION func(x,n) Function that outputs the n basis functions
! INTENT(INOUT)
! None
! INTENT(OUT)
! real(sp/dp), dimension(M) :: a fitted M parameters
! real(sp/dp) :: chisq Minimum chi^2
! real(sp/dp), dimension(M,M) :: v output for svdvar
! real(sp/dp), dimension(M) :: w output for svdvar
! INTENT(IN), OPTIONAL
! None
! INTENT(INOUT), OPTIONAL
! None
! INTENT(OUT), OPTIONAL
! None
! RESTRICTIONS
! None
! EXAMPLE
! ! fit polynomial with degree 5
! integer(i4b), parameter :: npt=100, npol=5
! integer(i4b) :: i
! real(dp) :: chisq
! real(dp), dimension(npol) :: a, w
! real(dp), dimension(npt) :: x, y, sig
! real(dp), dimension(npol,npol) :: cvm, v
! call svdfit(x, y, sig, a, v, w, chisq, fpoly_dp)
! call svdvar(v, w, cvm)
! do i=1, npol
! write(*,'(1x,f12.6,a,f10.6)') a(i), ' +-', sqrt(cvm(i,i))
! end do
! LITERATURE
! Press WH, Teukolsky SA, Vetterling WT, & Flannery BP - Numerical Recipes in Fortran 90 -
! The Art of Parallel Scientific Computing, 2nd Edition, Volume 2 of Fortran Numerical Recipes,
! Cambridge University Press, UK, 1996
! HISTORY
! Written, Matthias Cuntz, Mar 2011
INTERFACE svdfit
MODULE PROCEDURE svdfit_sp, svdfit_dp
END INTERFACE svdfit
! ------------------------------------------------------------------
! NAME
! svdvar
! PURPOSE
! Calculates the variance/covariance for the singular value decomposition estimates.
! From the numerical recipes documentation:
! To evaluate the covariance matrix cvm of the fit for ma parameters obtained by svdfit,
! call this routine with matrices v,w as returned from svdfit. np,ncvm give the physical
! dimensions of v,w,cvm as indicated.
! CALLING SEQUENCE
! call svdvar(v, w, cvm)
! INTENT(IN)
! real(sp/dp), dimension(M,M) :: v output from svdfit
! real(sp/dp), dimension(M) :: w output from svdfit
! INTENT(INOUT)
! None
! INTENT(OUT)
! real(sp/dp), dimension(M,M) :: cvm covariance matrix for M fitted parameters
! INTENT(IN), OPTIONAL
! None
! INTENT(INOUT), OPTIONAL
! None
! INTENT(OUT), OPTIONAL
! None
! RESTRICTIONS
! None
! EXAMPLE
! ! fit polynomial with degree 5
! integer(i4b), parameter :: npt=100, npol=5
! integer(i4b) :: i
! real(dp) :: chisq
! real(dp), dimension(npol) :: a, w
! real(dp), dimension(npt) :: x, y, sig
! real(dp), dimension(npol,npol) :: cvm, v
! call svdfit(x, y, sig, a, v, w, chisq, fpoly_dp)
! call svdvar(v, w, cvm)
! do i=1, npol
! write(*,'(1x,f12.6,a,f10.6)') a(i), ' +-', sqrt(cvm(i,i))
! end do
! LITERATURE
! Press WH, Teukolsky SA, Vetterling WT, & Flannery BP - Numerical Recipes in Fortran 90 -
! The Art of Parallel Scientific Computing, 2nd Edition, Volume 2 of Fortran Numerical Recipes,
! Cambridge University Press, UK, 1996
! HISTORY
! Written, Matthias Cuntz, Mar 2011
INTERFACE svdvar
MODULE PROCEDURE svdvar_sp, svdvar_dp
END INTERFACE svdvar
! ------------------------------------------------------------------
PRIVATE
! ------------------------------------------------------------------
! Private interfaces, mostly from numerical recipes
INTERFACE geop
MODULE PROCEDURE geop_dp, geop_i4, geop_sp, geop_v_dp, geop_v_sp
END INTERFACE geop
INTERFACE outerprod
MODULE PROCEDURE outerprod_dp, outerprod_sp
END INTERFACE outerprod
INTERFACE pythag
MODULE PROCEDURE pythag_dp, pythag_sp
END INTERFACE pythag
INTERFACE svdksb
MODULE PROCEDURE svdksb_dp, svdksb_sp
END INTERFACE svdksb
INTERFACE svdcmp
MODULE PROCEDURE svdcmp_dp, svdcmp_sp
END INTERFACE svdcmp
INTERFACE vabs
MODULE PROCEDURE vabs_dp, vabs_sp
END INTERFACE vabs
! for geop
INTEGER(i4), PARAMETER :: NPAR_GEOP=4, NPAR2_GEOP=2
! ------------------------------------------------------------------
CONTAINS
! ------------------------------------------------------------------
SUBROUTINE fitfun_dp(xin, yin, sigin, a, func, mask, fit, chisq, siga)
IMPLICIT NONE
REAL(dp), DIMENSION(:), INTENT(IN) :: xin
REAL(dp), DIMENSION(:), INTENT(IN) :: yin
REAL(dp), DIMENSION(:), INTENT(IN) :: sigin
REAL(dp), DIMENSION(:), INTENT(OUT) :: a
INTERFACE
FUNCTION func(x,n)
USE mo_kind, ONLY: i4, dp
IMPLICIT NONE
REAL(dp), INTENT(IN) :: x
INTEGER(i4), INTENT(IN) :: n
REAL(dp), DIMENSION(n) :: func
END FUNCTION func
END INTERFACE
LOGICAL, DIMENSION(:), OPTIONAL, INTENT(IN) :: mask
REAL(dp), DIMENSION(size(xin)), OPTIONAL, INTENT(OUT) :: fit
REAL(dp), OPTIONAL, INTENT(OUT) :: chisq
REAL(dp), DIMENSION(size(a)), OPTIONAL, INTENT(OUT) :: siga
LOGICAL, DIMENSION(size(xin)) :: maske
REAL(dp), DIMENSION(size(a)) :: w
REAL(dp), DIMENSION(size(a),size(a)) :: cvm, v
REAL(dp), DIMENSION(:), ALLOCATABLE :: x, y, sig
REAL(dp), DIMENSION(:), ALLOCATABLE :: out
REAL(dp) :: ichisq
INTEGER :: i, nmask, na
if (size(xin) /= size(yin)) stop 'Error fitfun_dp: size(x) /= size(y)'
if (size(sigin) /= size(yin)) stop 'Error fitfun_dp: size(sig) /= size(y)'
maske(:) = .true.
if (present(mask)) then
if (size(mask) /= size(xin)) stop 'Error fitfun_dp: size(mask) /= size(x)'
maske = mask
endif
nmask = count(maske)
allocate(x(nmask))
allocate(y(nmask))
allocate(sig(nmask))
x(1:nmask) = pack(xin,maske)
y(1:nmask) = pack(yin,maske)
sig(1:nmask) = pack(sigin,maske)
call svdfit(x, y, sig, a, v, w, ichisq, func)
if (present(chisq)) chisq = ichisq
na = size(a)
if (present(siga)) then
call svdvar(v, w, cvm)
do i=1, na
siga(i) = sqrt(cvm(i,i))
end do
end if
if (present(fit)) then
allocate(out(nmask))
do i=1, nmask
out(i) = dot_product(a, real(func(real(x(i),dp),na),dp))
end do
fit = unpack(out,maske,0.0_dp)
deallocate(out)
endif
deallocate(x)
deallocate(y)
deallocate(sig)
END SUBROUTINE fitfun_dp
SUBROUTINE fitfun_sp(xin, yin, sigin, a, func, mask, fit, chisq, siga)
IMPLICIT NONE
REAL(sp), DIMENSION(:), INTENT(IN) :: xin
REAL(sp), DIMENSION(:), INTENT(IN) :: yin
REAL(sp), DIMENSION(:), INTENT(IN) :: sigin
REAL(sp), DIMENSION(:), INTENT(OUT) :: a
INTERFACE
FUNCTION func(x,n)
USE mo_kind, ONLY: i4, sp
IMPLICIT NONE
REAL(sp), INTENT(IN) :: x
INTEGER(i4), INTENT(IN) :: n
REAL(sp), DIMENSION(n) :: func
END FUNCTION func
END INTERFACE
LOGICAL, DIMENSION(:), OPTIONAL, INTENT(IN) :: mask
REAL(sp), DIMENSION(size(xin)), OPTIONAL, INTENT(OUT) :: fit
REAL(sp), OPTIONAL, INTENT(OUT) :: chisq
REAL(sp), DIMENSION(size(a)), OPTIONAL, INTENT(OUT) :: siga
LOGICAL, DIMENSION(size(xin)) :: maske
REAL(sp), DIMENSION(size(a)) :: w
REAL(sp), DIMENSION(size(a),size(a)) :: cvm, v
REAL(sp), DIMENSION(:), ALLOCATABLE :: x, y, sig
REAL(sp), DIMENSION(:), ALLOCATABLE :: out
REAL(sp) :: ichisq
INTEGER :: i, nmask, na
if (size(xin) /= size(yin)) stop 'Error fitfun_sp: size(x) /= size(y)'
if (size(sigin) /= size(yin)) stop 'Error fitfun_sp: size(sig) /= size(y)'
maske(:) = .true.
if (present(mask)) then
if (size(mask) /= size(xin)) stop 'Error fitfun_sp: size(mask) /= size(x)'
maske = mask
endif
nmask = count(maske)
allocate(x(nmask))
allocate(y(nmask))
allocate(sig(nmask))
x(1:nmask) = pack(xin,maske)
y(1:nmask) = pack(yin,maske)
sig(1:nmask) = pack(sigin,maske)
call svdfit(x, y, sig, a, v, w, ichisq, func)
if (present(chisq)) chisq = ichisq
na = size(a)
if (present(siga)) then
call svdvar(v, w, cvm)
do i=1, na
siga(i) = sqrt(cvm(i,i))
end do
end if
if (present(fit)) then
allocate(out(nmask))
do i=1, nmask
out(i) = dot_product(a, real(func(real(x(i),sp),na),sp))
end do
fit = unpack(out,maske,0.0_sp)
deallocate(out)
endif
deallocate(x)
deallocate(y)
deallocate(sig)
END SUBROUTINE fitfun_sp
! ------------------------------------------------------------------
FUNCTION fpoly_dp(x,n)
IMPLICIT NONE
REAL(dp), INTENT(IN) :: x
INTEGER(i4), INTENT(IN) :: n
REAL(dp), DIMENSION(n) :: fpoly_dp
fpoly_dp = geop(1.0_dp,x,n)
END FUNCTION fpoly_dp
FUNCTION fpoly_sp(x,n)
IMPLICIT NONE
REAL(sp), INTENT(IN) :: x
INTEGER(i4), INTENT(IN) :: n
REAL(sp), DIMENSION(n) :: fpoly_sp
fpoly_sp = geop(1.0_sp,x,n)
END FUNCTION fpoly_sp
! ------------------------------------------------------------------
! From numerical recipes documentation
! Returns an array of length n containing a geometric progression whose first
! value is first and whose multiplier is factor. If first and factor
! have rank greater than zero, returns an array of one larger rank, with the
! last subscript having size n and indexing the progression.
FUNCTION geop_dp(first,factor,n)
REAL(dp), INTENT(IN) :: first,factor
INTEGER(i4), INTENT(IN) :: n
REAL(dp), DIMENSION(n) :: geop_dp
INTEGER(i4) :: k,k2
REAL(dp) :: temp
if (n > 0) geop_dp(1)=first
if (n <= NPAR_GEOP) then
do k=2,n
geop_dp(k)=geop_dp(k-1)*factor
end do
else
do k=2,NPAR2_GEOP
geop_dp(k)=geop_dp(k-1)*factor
end do
temp=factor**NPAR2_GEOP
k=NPAR2_GEOP
do
if (k >= n) exit
k2=k+k
geop_dp(k+1:min(k2,n))=temp*geop_dp(1:min(k,n-k))
temp=temp*temp
k=k2
end do
end if
END FUNCTION geop_dp
FUNCTION geop_i4(first,factor,n)
INTEGER(i4), INTENT(IN) :: first, factor, n
INTEGER(i4), DIMENSION(n) :: geop_i4
INTEGER(i4) :: k, k2, temp
if (n > 0) geop_i4(1)=first
if (n <= NPAR_GEOP) then
do k=2,n
geop_i4(k)=geop_i4(k-1)*factor
end do
else
do k=2,NPAR2_GEOP
geop_i4(k)=geop_i4(k-1)*factor
end do
temp=factor**NPAR2_GEOP
k=NPAR2_GEOP
do
if (k >= n) exit
k2=k+k
geop_i4(k+1:min(k2,n))=temp*geop_i4(1:min(k,n-k))
temp=temp*temp
k=k2
end do
end if
END FUNCTION geop_i4
FUNCTION geop_sp(first,factor,n)
REAL(sp), INTENT(IN) :: first,factor
INTEGER(i4), INTENT(IN) :: n
REAL(sp), DIMENSION(n) :: geop_sp
INTEGER(i4) :: k,k2
REAL(sp) :: temp
if (n > 0) geop_sp(1)=first
if (n <= NPAR_GEOP) then
do k=2,n
geop_sp(k)=geop_sp(k-1)*factor
end do
else
do k=2,NPAR2_GEOP
geop_sp(k)=geop_sp(k-1)*factor
end do
temp=factor**NPAR2_GEOP
k=NPAR2_GEOP
do
if (k >= n) exit
k2=k+k
geop_sp(k+1:min(k2,n))=temp*geop_sp(1:min(k,n-k))
temp=temp*temp
k=k2
end do
end if
END FUNCTION geop_sp
FUNCTION geop_v_dp(first,factor,n)
REAL(dp), DIMENSION(:), INTENT(IN) :: first, factor
INTEGER(i4), INTENT(IN) :: n
REAL(dp), DIMENSION(size(first),n) :: geop_v_dp
INTEGER(i4) :: k,k2
REAL(dp), DIMENSION(size(first)) :: temp
if (n > 0) geop_v_dp(:,1)=first(:)
if (n <= NPAR_GEOP) then
do k=2,n
geop_v_dp(:,k)=geop_v_dp(:,k-1)*factor(:)
end do
else
do k=2,NPAR2_GEOP
geop_v_dp(:,k)=geop_v_dp(:,k-1)*factor(:)
end do
temp=factor**NPAR2_GEOP
k=NPAR2_GEOP
do
if (k >= n) exit
k2=k+k
geop_v_dp(:,k+1:min(k2,n))=geop_v_dp(:,1:min(k,n-k))*&
spread(temp,2,size(geop_v_dp(:,1:min(k,n-k)),2))
temp=temp*temp
k=k2
end do
end if
END FUNCTION geop_v_dp
FUNCTION geop_v_sp(first,factor,n)
REAL(sp), DIMENSION(:), INTENT(IN) :: first, factor
INTEGER(i4), INTENT(IN) :: n
REAL(sp), DIMENSION(size(first),n) :: geop_v_sp
INTEGER(i4) :: k,k2
REAL(sp), DIMENSION(size(first)) :: temp
if (n > 0) geop_v_sp(:,1)=first(:)
if (n <= NPAR_GEOP) then
do k=2,n
geop_v_sp(:,k)=geop_v_sp(:,k-1)*factor(:)
end do
else
do k=2,NPAR2_GEOP
geop_v_sp(:,k)=geop_v_sp(:,k-1)*factor(:)
end do
temp=factor**NPAR2_GEOP
k=NPAR2_GEOP
do
if (k >= n) exit
k2=k+k
geop_v_sp(:,k+1:min(k2,n))=geop_v_sp(:,1:min(k,n-k))*&
spread(temp,2,size(geop_v_sp(:,1:min(k,n-k)),2))
temp=temp*temp
k=k2
end do
end if
END FUNCTION geop_v_sp
! ------------------------------------------------------------------
! From numerical recipes documentation
! Returns a matrix that is the outer product of two vectors.
! a*b with size(a)/=size(b)
FUNCTION outerprod_dp(a,b)
REAL(dp), DIMENSION(:), INTENT(IN) :: a, b
REAL(dp), DIMENSION(size(a),size(b)) :: outerprod_dp
outerprod_dp = spread(a,dim=2,ncopies=size(b)) * &
spread(b,dim=1,ncopies=size(a))
END FUNCTION outerprod_dp
FUNCTION outerprod_sp(a,b)
REAL(sp), DIMENSION(:), INTENT(IN) :: a, b
REAL(sp), DIMENSION(size(a),size(b)) :: outerprod_sp
outerprod_sp = spread(a,dim=2,ncopies=size(b)) * &
spread(b,dim=1,ncopies=size(a))
END FUNCTION outerprod_sp
! ------------------------------------------------------------------
! sqrt(a**2+b**2)
FUNCTION pythag_dp(a,b)
IMPLICIT NONE
REAL(dp), INTENT(IN) :: a,b
REAL(dp) :: pythag_dp
REAL(dp) :: absa, absb
absa = abs(a)
absb = abs(b)
if (absa > absb) then
pythag_dp = absa*sqrt(1.0_dp+(absb/absa)**2)
else
if (absb < tiny(1.0_dp)) then
pythag_dp = 0.0_dp
else
pythag_dp = absb*sqrt(1.0_dp+(absa/absb)**2)
end if
end if
END FUNCTION pythag_dp
FUNCTION pythag_sp(a,b)
IMPLICIT NONE
REAL(sp), INTENT(IN) :: a,b
REAL(sp) :: pythag_sp
REAL(sp) :: absa, absb
absa = abs(a)
absb = abs(b)
if (absa > absb) then
pythag_sp = absa*sqrt(1.0_sp+(absb/absa)**2)
else
if (absb < tiny(1.0_sp)) then
pythag_sp = 0.0_sp
else
pythag_sp = absb*sqrt(1.0_sp+(absa/absb)**2)
end if
end if
END FUNCTION pythag_sp
! ------------------------------------------------------------------
FUNCTION polyfit_dp(vxin, vyin, d, mask)
IMPLICIT NONE
real(dp), dimension(:), intent(in) :: vxin, vyin
integer, intent(in) :: d
LOGICAL, DIMENSION(:), OPTIONAL, INTENT(IN) :: mask
real(dp), dimension(d+1) :: polyfit_dp
real(dp), dimension(:,:), allocatable :: X
real(dp), dimension(:,:), allocatable :: XT
real(dp), dimension(:,:), allocatable :: XTX
integer :: i
integer :: n, lda, lwork
integer :: info
integer, dimension(:), allocatable :: ipiv
real(dp), dimension(:), allocatable :: work
real(dp), dimension(:), allocatable :: vx, vy
LOGICAL, DIMENSION(size(vxin)) :: maske
integer :: nmask
EXTERNAL :: DGETRF, DGETRI
maske(:) = .true.
if (present(mask)) then
if (size(mask) /= size(vxin)) stop 'Error polyfit_dp: size(mask) /= size(vxin)'
maske = mask
endif
nmask = count(maske)
allocate(vx(nmask))
allocate(vy(nmask))
vx(1:nmask) = pack(vxin,maske)
vy(1:nmask) = pack(vyin,maske)
n = d+1
lda = n
lwork = n
allocate(ipiv(n))
allocate(work(lwork))
allocate(XT(n,size(vx)))
allocate(X(size(vx),n))
allocate(XTX(n,n))
! prepare the matrix
X(1:nmask,1) = 1.0_dp
if (n > 1) X(1:nmask,2) = vx
do i=3, n
X(1:nmask,i) = vx(1:nmask)**(i-1)
end do
XT = transpose(X)
XTX = matmul(XT,X)
! calls to LAPACK subs DGETRF and DGETRI
call DGETRF(n, n, XTX, lda, ipiv, info)
if ( info /= 0 ) stop 'Problem polyfit_dp 1'
call DGETRI(n, XTX, lda, ipiv, work, lwork, info)
if ( info /= 0 ) stop 'Problem polyfit_dp 2'
polyfit_dp = matmul(matmul(XTX,XT), vy)
deallocate(ipiv)
deallocate(work)
deallocate(X)
deallocate(XT)