-
Notifications
You must be signed in to change notification settings - Fork 10
/
Copy pathminimizer_engine.f90
1515 lines (1099 loc) · 45.5 KB
/
minimizer_engine.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!
! Copyright 2011 Sebastian Heimann
!
! Licensed under the Apache License, Version 2.0 (the "License");
! you may not use this file except in compliance with the License.
! You may obtain a copy of the License at
!
! http://www.apache.org/licenses/LICENSE-2.0
!
! Unless required by applicable law or agreed to in writing, software
! distributed under the License is distributed on an "AS IS" BASIS,
! WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
! See the License for the specific language governing permissions and
! limitations under the License.
!
module minimizer_engine
! this thing represents a source - greensfunction - receiver - seismogram - misfit setup
use constants
use util
use unit
use orthodrome
use gfdb
use source
use seismogram
use seismogram_io
use sparse_trace
use receiver
use better_varying_string
use read_table
use comparator
use piecewise_linear_function
use omp_lib
implicit none
private
public set_database
public set_local_interpolation
public set_spacial_undersampling
public set_ref_seismograms
public set_source_location
public set_source_constraints
public set_source_crustal_thickness_limit
public get_source_crustal_thickness
public set_source_params
public set_source_params_mask
public set_source_subparams
public set_source_subparams_limits
public set_effective_dt
public set_receivers
public switch_receiver
public set_misfit_method
public set_misfit_filter
public set_misfit_taper
public set_synthetics_factor
public minimize_lm
public output_seismograms
public output_seismogram_spectra
public output_source_model
public cleanup_minimizer
public get_source_subparams
public get_distances
public get_global_misfit
public get_misfits
public get_floating_shifts
public get_peak_amplitudes
public get_arias_intensities
public get_principal_axes
public output_cross_correlations
public shift_ref_seismogram
public autoshift_ref_seismogram
public get_cached_traces_memory
public set_cached_traces_memory_limit
public set_floating_shiftrange
type(t_psm), save :: psm
real :: effective_dt = 1.
type(t_tdsm) :: tdsm
type(t_receiver), dimension(:), allocatable :: receivers
type(t_gfdb), save :: db
integer :: misfit_method = L2NORM
real :: misfit
logical :: interpolate = .false.
integer :: xundersample = 1
integer :: zundersample = 1
real, dimension(:), allocatable :: g_subparam_mins
real, dimension(:), allocatable :: g_subparam_maxs
logical :: database_inited = .false.
logical :: ref_probes_inited = .false.
logical :: syn_probes_inited = .false.
logical :: receivers_inited = .false.
logical :: source_inited = .false.
logical :: source_location_inited = .false.
logical :: seismograms_inited = .false.
logical :: misfits_inited = .false.
logical :: database_dirty = .true.
logical :: ref_probes_dirty = .true.
logical :: syn_probes_dirty = .true.
logical :: receivers_dirty = .true.
logical :: source_dirty = .true.
logical :: source_location_dirty = .true.
logical :: seismograms_dirty = .true.
logical :: misfits_dirty = .true.
integer :: iterations ! needed by lm_forward_step and minmize_lm
contains
subroutine set_database( db_path, nipx, nipz, ok )
type(varying_string), intent(in) :: db_path
integer, intent(in) :: nipx, nipz
logical, intent(out) :: ok
real :: olddt
ok = .true.
olddt = db%dt
call gfdb_init(db, db_path, nipx=nipx, nipz=nipz ) ! dies program if not successful
if (database_inited) then
if (db%dt /= olddt) then
call warn("sampling rate of new database is different than that of the old database")
call warn("deinitializing receivers...")
call cleanup_receivers()
end if
end if
database_inited = .true.
call dirtyfy_database()
end subroutine
subroutine set_local_interpolation( state )
logical, intent(in) :: state
interpolate = state
call dirtyfy_seismograms()
end subroutine
subroutine set_spacial_undersampling( xunder, zunder, ok )
integer, intent(in) :: xunder, zunder
logical, intent(out) :: ok
ok = .true.
if (xunder < 1 .or. zunder < 1) then
call error( "invalid undersampling value" )
ok = .false.
end if
xundersample = xunder
zundersample = zunder
call dirtyfy_database()
end subroutine
subroutine set_receivers( receiversfn, has_depth, answer, ok )
! load receiver list from file
type(varying_string), intent(in) :: receiversfn
logical, intent(in) :: has_depth
type(varying_string), intent(out) :: answer
logical, intent(out) :: ok
character(len=1024) :: components, str
integer :: ireceiver, nreceivers
integer :: iostat, nskip, iunit
integer :: nwords
type(t_geo_coords) :: origin
real :: depth
! this subroutine will be a pain in fortran
answer = ''
ok = .true.
call cleanup_receivers()
! open the file
call claim_unit( iunit )
open( unit=iunit, file=char(receiversfn), status="old", action="read", iostat=iostat )
if (iostat /= 0) then
call release_unit(iunit)
call error( "can't open file " // receiversfn )
ok = .false.
return
end if
! determine how many receivers have to be allocated and allocate them
nreceivers = count_non_comment_lines( iunit, iostat )
if (iostat /= 0) then
close(iunit)
call release_unit(iunit)
call error( "error occured while counting non comment lines in " // receiversfn )
ok = .false.
return
end if
allocate( receivers( nreceivers ) )
! now really read the file and initialize the receivers
ireceiver = 0
line_loop : do
nskip = skip_comments(iunit,iostat)
if (iostat == IOSTAT_EOF .or. nskip < 0) exit line_loop
if (iostat /= 0) then
ok = .false.
call error("io error occured while skipping comments on file " // receiversfn )
exit line_loop
end if
read (iunit,"(a)",iostat=iostat) str
if (iostat == IOSTAT_EOF) exit line_loop
if (iostat /= 0) then
ok = .false.
call error("io error occured while reading " // receiversfn // " at receiver no " // (ireceiver+1) )
exit line_loop
end if
components = ""
depth = 0.0
nwords = count_words(str)
if (has_depth .and. nwords == 4) then
read (str,*,iostat=iostat) origin%lat, origin%lon, depth, components
else if (has_depth .and. nwords == 3) then
read (str,*,iostat=iostat) origin%lat, origin%lon, depth
else if (.not. has_depth .and. nwords == 3) then
read (str,*,iostat=iostat) origin%lat, origin%lon, components
else if (.not. has_depth .and. nwords == 2) then
read (str,*,iostat=iostat) origin%lat, origin%lon
else
ok = .false.
call error("expected two or three words at receiver no " // (ireceiver+1) // &
" while reading " // receiversfn)
exit line_loop
end if
if (iostat /= 0) then
ok = .false.
call error("io error occured while reading " // receiversfn // " at receiver no " // (ireceiver+1) )
exit line_loop
end if
ireceiver = ireceiver + 1
if (ireceiver > nreceivers) exit line_loop
call receiver_init(receivers(ireceiver), d2r(origin), depth, components, db%dt, ok )
call receiver_set_ids(receivers(ireceiver), var_str(''), var_str(ireceiver), var_str(''))
if (.not. ok) then
call error("initializing receiver failed: possibly a forbidden "// &
"combination of receiver components has been given at receiver no. "// ireceiver)
exit line_loop
end if
end do line_loop
close(iunit)
call release_unit(iunit)
! cleanup and say goodbye if something went wrong
if (.not. ok) then
call cleanup_receivers()
return
end if
receivers_inited = .true.
ref_probes_inited = .false.
call dirtyfy_receivers()
end subroutine
subroutine switch_receiver( ireceiver, newstate, ok )
integer, intent(in) :: ireceiver
logical, intent(in) :: newstate
logical, intent(out) :: ok
integer :: nreceivers
call update_receivers( ok )
if (.not. ok) return
nreceivers = size(receivers)
if (ireceiver < 1 .or. nreceivers < ireceiver) then
ok = .false.
call error( 'receiver index out of range' )
return
end if
call receiver_set_enabled( receivers(ireceiver), newstate )
call dirtyfy_receivers()
end subroutine
subroutine set_ref_seismograms( reffnbase, refformat, ok )
type(varying_string), intent(in) :: reffnbase, refformat
logical, intent(out) :: ok
! read a set of reference seismograms from ascii or sac files
integer :: ireceiver, nreceivers
type(varying_string) :: reffn
call update_database( ok )
if (.not. ok) then
call warn("have to set database prior to setting ref seismograms")
return
end if
call update_source_location( ok )
if (.not. ok) then
call warn("have to set source location prior to setting ref seismograms")
return
end if
call update_receivers( ok )
if (.not. ok) return
nreceivers = size(receivers)
do ireceiver=1,nreceivers
reffn = reffnbase // "-" // ireceiver
call receiver_set_ref_seismogram( receivers(ireceiver), reffn, refformat, psm%ref_time, ok )
if (.not. ok) then
ref_probes_inited = .false.
return
end if
end do
ref_probes_inited = .true.
call dirtyfy_ref_probes()
end subroutine
subroutine shift_ref_seismogram( ireceiver, shift, ok )
integer, intent(in) :: ireceiver
real, intent(in) :: shift
logical, intent(out) :: ok
integer :: nreceivers
integer :: ishift
call update_ref_probes( ok )
if (.not. ok) return
nreceivers = size(receivers)
if (ireceiver < 1 .or. nreceivers < ireceiver) then
ok = .false.
call error( 'receiver index out of range' )
return
end if
ishift = int(nint(shift/db%dt))
call receiver_shift_ref_seismogram( receivers(ireceiver), ishift )
call dirtyfy_ref_probes()
end subroutine
subroutine autoshift_ref_seismogram( ireceiver, shiftrange, shifts, ok )
integer, intent(in) :: ireceiver
real, dimension(2), intent(in) :: shiftrange
real, dimension(:), allocatable, intent(inout) :: shifts
logical, intent(out) :: ok
integer :: nreceivers, irec, ishift
integer, dimension(2) :: ishiftrange
ok = .true.
if (allocated(shifts)) deallocate(shifts)
call update_misfits(ok)
if (.not. ok) return
nreceivers = size(receivers)
ishiftrange(:) = int(nint(shiftrange(:)/db%dt))
if (ireceiver == 0) then ! apply to all
allocate( shifts(nreceivers) )
do irec=1,nreceivers
call receiver_autoshift_ref_seismogram( receivers(irec), ishiftrange, ishift )
shifts(irec) = ishift*db%dt
end do
else
if (ireceiver < 1 .or. nreceivers < ireceiver) then
ok = .false.
call error( 'receiver index out of range' )
return
end if
allocate(shifts(1))
call receiver_autoshift_ref_seismogram( receivers(ireceiver), ishiftrange, ishift )
shifts(1) = ishift*db%dt
end if
call dirtyfy_ref_probes()
end subroutine
subroutine set_floating_shiftrange( ireceiver, shiftrange, ok )
integer, intent(in) :: ireceiver
real, dimension(2), intent(in) :: shiftrange
logical, intent(out) :: ok
integer :: nreceivers, irec
integer, dimension(2) :: ishiftrange
ok = .true.
nreceivers = size(receivers)
ishiftrange(:) = int(nint(shiftrange(:)/db%dt))
if (ireceiver == 0) then ! apply to all
do irec=1,nreceivers
receivers(irec)%floating_shiftrange = ishiftrange
end do
else
if (ireceiver < 1 .or. nreceivers < ireceiver) then
ok = .false.
call error( 'receiver index out of range' )
return
end if
receivers(ireceiver)%floating_shiftrange = ishiftrange
end if
call dirtyfy_ref_probes()
end subroutine
subroutine set_source_location( lat, lon, ref_time )
real, intent(in) :: lat, lon
double precision, intent(in) :: ref_time
type(t_geo_coords) :: origin
origin%lat = lat
origin%lon = lon
call psm_set_origin_and_time(psm, origin, ref_time)
source_location_inited = .true.
call dirtyfy_source_location()
end subroutine
subroutine set_source_constraints( points, normals )
real, dimension(:,:), intent(in) :: points
real, dimension(:,:), intent(in) :: normals
call psm_set_constraints( psm, points, normals )
call dirtyfy_source()
end subroutine
subroutine set_source_crustal_thickness_limit( thickness_limit )
real, intent(in) :: thickness_limit
call psm_set_crustal_thickness_limit(psm, thickness_limit)
call dirtyfy_source()
end subroutine
subroutine get_source_crustal_thickness( thickness, ok )
real, intent(out) :: thickness
logical, intent(out) :: ok
call update_source_location( ok )
if (.not. ok) return
call psm_get_crustal_thickness(psm, thickness)
end subroutine
subroutine set_source_params( sourcetype, sourceparams, ok )
integer, intent(in) :: sourcetype
real, dimension(:),intent(in) :: sourceparams
logical, intent(out) :: ok
logical :: only_moment_changed
call update_source_location( ok )
if (.not. ok) return
if (source_inited .and. sourcetype == psm%sourcetype) then
if (all(psm%params == sourceparams)) return ! this source has already been set
end if
call psm_set(psm, sourcetype, sourceparams, only_moment_changed_=only_moment_changed )
if (only_moment_changed) then
call dirtyfy_syn_probes()
else
source_inited = .true.
call dirtyfy_source()
end if
end subroutine
subroutine set_source_params_mask( params_mask , ok)
logical, dimension(:), intent(in) :: params_mask
logical, intent(out) :: ok
call update_source( ok )
if (.not. ok) return
if (size(psm%params) /= size(params_mask)) then
ok = .false.
call error( 'wrong number of elements in mask' )
return
end if
call psm_set_mask( psm, params_mask )
call reset_subparam_limits()
end subroutine
subroutine set_source_subparams( subparams, ok )
real, dimension(:), intent(in) :: subparams
logical, intent(out) :: ok
ok = .true.
if (.not. source_inited) then
ok = .false.
call error( 'source parameters must be set prior to setting parameter subset' )
return
end if
if (count(psm%params_mask) /= size(subparams)) then
ok = .false.
call error( 'wrong number of subparams' )
return
end if
call psm_set_subparams( psm, subparams )
call dirtyfy_source()
end subroutine
subroutine reset_subparam_limits()
if (allocated(g_subparam_mins)) then
deallocate(g_subparam_mins)
end if
if (allocated(g_subparam_maxs)) then
deallocate(g_subparam_maxs)
end if
end subroutine
subroutine set_source_subparams_limits( subparam_mins, subparam_maxs, ok )
real, dimension(:), intent(in) :: subparam_mins
real, dimension(:), intent(in) :: subparam_maxs
logical, intent(out) :: ok
integer :: nsubparams
nsubparams = count(psm%params_mask)
if (nsubparams /= size(subparam_mins)) then
ok = .false.
call error( 'wrong number of subparam_mins' )
return
end if
if (nsubparams /= size(subparam_maxs)) then
ok = .false.
call error( 'wrong number of subparam_maxs' )
return
end if
call reset_subparam_limits()
allocate(g_subparam_mins(nsubparams))
allocate(g_subparam_maxs(nsubparams))
g_subparam_mins(:) = subparam_mins(:)
g_subparam_maxs(:) = subparam_maxs(:)
end subroutine
subroutine set_effective_dt( effective_dt_ )
real, intent(in) :: effective_dt_
effective_dt = effective_dt_
call dirtyfy_source()
end subroutine
subroutine set_misfit_method( norm )
integer, intent(in) :: norm
misfit_method = norm
call dirtyfy_misfits()
end subroutine
subroutine set_misfit_filter( ireceiver, x, y, ok )
! set plf filter on one or all receivers
integer, intent(in) :: ireceiver
real, dimension(:), intent(in) :: x,y
logical, intent(out) :: ok
integer :: irec, nreceivers
type(t_plf) :: filter
call update_receivers( ok )
if (.not. ok) return
nreceivers = size(receivers)
if (ireceiver < 0 .or. nreceivers < ireceiver) then
ok = .false.
call error( 'receiver index out of range' )
return
end if
call plf_make( filter, x, y )
if (ireceiver .eq. 0) then
do irec=1,nreceivers
call receiver_set_filter( receivers(irec), filter )
end do
else
call receiver_set_filter( receivers(ireceiver), filter )
end if
call plf_destroy( filter )
call dirtyfy_misfits()
end subroutine
subroutine set_misfit_taper( ireceiver, x, y, ok )
integer, intent(in) :: ireceiver
real, dimension(:), intent(in) :: x,y
logical, intent(out) :: ok
integer :: nreceivers
type(t_plf) :: taper
!call update_ref_probes( ok )
!if (.not. ok) return
call update_receivers( ok )
if (.not. ok) return
nreceivers = size(receivers)
if (ireceiver < 1 .or. nreceivers < ireceiver) then
ok = .false.
call error( 'receiver index out of range' )
return
end if
call plf_make( taper, x, y )
call receiver_set_taper( receivers(ireceiver), taper )
call plf_destroy( taper )
call dirtyfy_misfits()
end subroutine
subroutine set_synthetics_factor( factor )
real, intent(in) :: factor
integer :: i
do i=1,size(receivers)
call receiver_set_synthetics_factor( receivers(i), factor )
end do
call dirtyfy_misfits()
end subroutine
pure function get_nmisfits(receivers)
integer :: get_nmisfits
type(t_receiver), dimension(:), intent(in) :: receivers
integer :: i
get_nmisfits = 0
do i=1,size(receivers)
get_nmisfits = get_nmisfits + receivers(i)%ncomponents
end do
end function
subroutine minimize_lm( info, iterations_, misfit_, ok )
integer, intent(out) :: info, iterations_
real, intent(out) :: misfit_
logical, intent(out) :: ok
call update_misfits( ok )
if (.not. ok) return
call minimize_lm_( info, iterations_, misfit_, ok )
end subroutine
subroutine minimize_lm_( info, iterations_, misfit_, ok )
integer, intent(out) :: info, iterations_
real, intent(out) :: misfit_
logical, intent(out) :: ok
real, dimension(:), allocatable :: subparams
real, dimension(get_nmisfits(receivers)) :: misfits
integer, dimension(count(psm%params_mask)) :: iwa
integer, dimension(get_nmisfits(receivers)*count(psm%params_mask)+5*count(psm%params_mask)+get_nmisfits(receivers)) :: wa
integer :: lwa, nsubparams
integer :: nmisfits
real :: tol
real,dimension(count(psm%params_mask)) :: diag
integer :: maxfev,mode,mp5n,nfev,nprint
real :: epsfcn,ftol,gtol,xtol,factor
real :: spmpar ! external function (sminpack)
ok = .true.
lwa = size(wa)
nmisfits = size(misfits)
call psm_get_subparams( psm, subparams, normalized_=.true. )
nsubparams = size(subparams)
! set tol to the square root of the machine precision.
! unless high precision solutions are required,
! this is the recommended setting.
tol = sqrt(spmpar(1))
iterations = 0
info = 0
if (nsubparams .le. 0 .or. nmisfits .lt. nsubparams .or. tol .lt. 0. &
.or. lwa .lt. nmisfits*nsubparams + 5*nsubparams + nmisfits) &
call die("something went wrong in minimize_lm")
maxfev = 500*(nsubparams + 1)
ftol = tol
xtol = tol
gtol = 0.
epsfcn = 0.
mode = 2
nprint = 0
factor = 0.01
mp5n = nmisfits + 5*nsubparams
diag(:) = 1.
call lmdif(lm_forward_step,nmisfits,nsubparams,subparams,misfits, &
ftol,xtol,gtol,maxfev,epsfcn,diag, &
mode,factor,nprint,info,nfev,wa(mp5n+1),nmisfits,iwa, &
wa(nsubparams+1),wa(2*nsubparams+1),wa(3*nsubparams+1), &
wa(4*nsubparams+1),wa(5*nsubparams+1))
if (info .eq. 8) info = 4
call resize( subparams, 1, 0 )
iterations_ = iterations
misfit_ = misfit
end subroutine
subroutine lm_forward_step( nmisfits, nsubparams, subparams, misfits, iflag )
! called by lmdif
integer :: nmisfits, nsubparams, iflag
real :: subparams(nsubparams), misfits(nmisfits)
integer :: imisfit
integer :: ireceiver, icomponent, nreceivers
logical :: ok
real, dimension(:), allocatable :: subparams_nn
real :: penalty
integer :: i
real, dimension(:), allocatable :: subparams_norm
penalty = 0.0
if (allocated(g_subparam_mins) .and. allocated(g_subparam_maxs)) then
call psm_get_subparams_norm( psm, subparams_norm )
do i=1,nsubparams
if (subparams(i)*subparams_norm(i) < g_subparam_mins(i)) then
penalty = penalty + abs(subparams(i)*subparams_norm(i) &
- g_subparam_mins(i)) &
/ abs(g_subparam_maxs(i) - g_subparam_mins(i))
subparams(i) = g_subparam_mins(i) / subparams_norm(i)
end if
if (subparams(i)*subparams_norm(i) > g_subparam_maxs(i)) then
penalty = penalty + abs(subparams(i)*subparams_norm(i) &
- g_subparam_maxs(i)) &
/ abs(g_subparam_maxs(i) - g_subparam_mins(i))
subparams(i) = g_subparam_maxs(i) / subparams_norm(i)
end if
end do
deallocate(subparams_norm)
end if
call psm_set_subparams( psm, subparams, normalized_=.true. )
call dirtyfy_source()
! print not normalize subparams
call psm_get_subparams( psm, subparams_nn, normalized_=.false. )
! write (stderr,*) subparams_nn
deallocate( subparams_nn )
call update_misfits( ok )
if ( .not. ok ) then
iflag = -2
return
end if
! gather misfits which are stored in receiver objects in misfits array
imisfit = 1
nreceivers = size(receivers)
do ireceiver=1,nreceivers
do icomponent=1,receivers(ireceiver)%ncomponents
misfits(imisfit) = receivers(ireceiver)%misfits(icomponent) * &
(1.0 + penalty)
imisfit = imisfit + 1
end do
end do
iterations = iterations+1
end subroutine
subroutine discretize_source( ok )
logical, intent(out) :: ok
call inform("discretizing source")
! convert to discrete source
call psm_to_tdsm( psm, tdsm, effective_dt, ok )
end subroutine
subroutine calculate_seismograms()
! Seismograms are calculated for the current configuration
integer :: ireceiver, nreceivers
call inform("calculating seismograms")
nreceivers = size(receivers)
!$omp parallel shared(tdsm, db, receivers) private(ireceiver)
!$omp do schedule(dynamic)
do ireceiver=1,nreceivers
if (receivers(ireceiver)%enabled) then
call make_seismogram( tdsm, receivers(ireceiver), db, interpolate, xundersample, zundersample )
end if
end do
!$omp end do
!$omp end parallel
seismograms_inited = .true.
end subroutine
subroutine scale_seismograms()
! Put seismogram data into probes and apply moment and rise-time
integer :: irec
call inform("scaling seismograms " // psm%moment)
do irec=1,size(receivers)
call receiver_scaled_seismograms_to_probes(receivers(irec), psm%risetime, psm%moment)
end do
end subroutine
subroutine calculate_misfits()
! calculate distance between reference and synthetic seismograms
integer :: ireceiver, nreceivers
real :: misfit_norm_factor
call inform("calculating misfits")
nreceivers = size(receivers)
! let receivers do the misfit calculation and gather the global misfit
misfit = 0.
misfit_norm_factor = 0.
do ireceiver=1,nreceivers
call receiver_calculate_misfits( receivers(ireceiver), misfit_method )
misfit = misfit + sum( receivers(ireceiver)%misfits**2 )
misfit_norm_factor = misfit_norm_factor + sum( receivers(ireceiver)%misfits_norm_factors**2 )
end do
misfit = sqrt(misfit) / sqrt(misfit_norm_factor)
misfits_inited = .true.
end subroutine
subroutine output_source_model( filenamebase,ok )
type(varying_string), intent(in) :: filenamebase
logical, intent(out) :: ok
type(varying_string) :: infofn
integer :: ofile, icentroid, iostat
call update_source( ok )
if (.not. ok) return
infofn = filenamebase // "-psm.info"
call psm_write_info_file( psm, infofn )
infofn = filenamebase // "-tdsm.info"
call tdsm_write_info_file( tdsm, infofn )
! dump discrete source centroids to file
call claim_unit( ofile )
infofn = filenamebase // "-dsm.table"
open( unit=ofile, file=char(infofn), status='unknown', iostat=iostat )
if (iostat /= 0) call die( "failed to open output file: " // infofn )
do icentroid=1,size(tdsm%centroids)
write (unit=ofile, fmt=*) &
tdsm%centroids(icentroid)%north, tdsm%centroids(icentroid)%east, &
tdsm%centroids(icentroid)%depth, tdsm%centroids(icentroid)%time, tdsm%centroids(icentroid)%m
end do
close( ofile )
call release_unit( ofile )
end subroutine
subroutine output_seismograms( filenamebase, fileformat, which_probe, which_processing, ok )
type(varying_string), intent(in) :: filenamebase, fileformat
integer, intent(in) :: which_probe, which_processing
logical, intent(out) :: ok
type(varying_string) :: outfn
integer :: ireceiver, nreceivers
if (which_probe == SYNTHETICS) then
call update_syn_probes( ok )
if (.not. ok) return
else if (which_probe == REFERENCES) then
call update_ref_probes( ok )
if (.not. ok) return
else
ok = .false.
call error("output_seismograms(): which_probe has invalid value")
return
end if