-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathpetscinv.F90
4549 lines (3748 loc) · 188 KB
/
petscinv.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
! -------------------------------------------------------------------------
!
! PETScinv:
!
! Solves a tomographic linear system in parallel with KSP. Allows
! for solution of rectangular systmes using a parallel implementaion
! of lsqr, or the normal equations using various direct solvers
!
! This is a complete Fortran 2008 rewrite of Lapo Boschi's original
! tomography toolbox. This is research code and we give no warranty
! nor do we guaranty fitness for a particular purpose.
!
! Compile with "make" and run with "./run_petscinv""
!
! Copyright (c) 2013 - 2015 Ludwig Auer, [email protected]
!
! -------------------------------------------------------------------------
!============================================================
!
! To specify all global parameters
!
module global_param
implicit none
public
!***********************************
! petsc-headers for module_mesh
! petscsys.h - base PETSc routines
! petscvec.h - vectors
! petscmat.h - matrices
! petscksp.h - Krylov subspace methods
! petscpc.h - preconditioners
#include <finclude/petscsys.h>
#include <finclude/petscvec.h>
#include <finclude/petscmat.h>
#include <finclude/petscis.h>
#include <finclude/petscao.h>
#include <finclude/petscpc.h>
#include <finclude/petscksp.h>
#include <finclude/petscviewer.h>
#include <finclude/petscdraw.h>
!***********************************
real(8), parameter :: pi = 3.1415926535898D0
real(8), parameter :: deg2rad = pi / 180.d0
real(8), parameter :: rad2deg = 180.d0 / pi
real(8), parameter :: r_earth = 6371.d0
real(8), parameter :: deg2km = 111.32 ! 1 deg in km at the equator (approx.)
real(8), parameter :: ref_eqinc = 5.0 ! reference grid, largest block size
integer, parameter :: sp = selected_real_kind(6, 37)
integer, parameter :: dp = selected_real_kind(15, 307)
integer, parameter :: qp = selected_real_kind(33, 4931)
integer, parameter :: verbosity = 1 ! 1=standard, 2=a lot of output
integer, parameter :: nlaym = 50 ! maximum number of layers
integer, parameter :: n0max = 50000 ! maximum number of blocks per layer
integer, parameter :: n1max = n0max*nlaym ! maximum number of blocks per param
integer, parameter :: npmax = 4 ! maximum number of physical params
integer, parameter :: nlmax = 288 ! maximum number of latitude zones
integer, parameter :: matmx = 3000 ! maximal number of matrices to read
integer, parameter :: dmpmx = 500 ! maximal number of damping loops
integer, parameter :: nprmx = 100000 ! maximal number per row for temporary arrays
integer, parameter :: datmx = 5*(10**6) ! maximal number of data per submatrix
integer, parameter :: nperd = 30 ! maximal number of nonz in damping row
PetscErrorCode ierr
PetscBool flg
PetscLogDouble memory
PetscInt irun
PetscInt imat
contains
!========================================================
!
! converts an integer to a string
!
character(len=20) function int2str(k)
integer, intent(in) :: k
write(int2str, *) k
int2str = adjustl(int2str)
end function int2str
!========================================================
end module global_param
!============================================================
!============================================================
!
! module related to reading options from command line
!
module module_options
use global_param
implicit none
type,public :: opts
PetscScalar eqinc
PetscScalar eqinc_ref
PetscInt npars
PetscInt nlays
PetscBool adapt
PetscBool gvarr
PetscChar(256) type
PetscChar(256) sched
PetscChar(256) grdin
PetscChar(256) inpar
PetscChar(256) param
PetscChar(256) synth
PetscChar(256) modarr
PetscChar(256) refmod
PetscChar(256) format
PetscChar(256) projid
character(5) :: parameter(4) ! maximum of 4 parameters
contains
procedure :: parse_options
procedure :: tokenize_param
end type opts
contains
!========================================================
!
! reads options from command line
!
subroutine parse_options(this)
implicit none
class(opts) :: this
! Initialize some default variables
this%eqinc_ref=ref_eqinc ! defined in the very top
this%synth=''
this%modarr=''
call PetscOptionsGetInt(PETSC_NULL_CHARACTER,&
'-number_of_layers',this%nlays,flg,ierr)
call PetscOptionsGetReal(PETSC_NULL_CHARACTER,&
'-equatorial_increment',this%eqinc,flg,ierr)
call PetscOptionsGetString(PETSC_NULL_CHARACTER,&
'-inversion_parameters',this%param,flg,ierr)
call PetscOptionsGetString(PETSC_NULL_CHARACTER,&
'-reference_model',this%refmod,flg,ierr)
call PetscOptionsHasName(PETSC_NULL_CHARACTER,&
'-adaptive_grid',this%adapt,ierr)
call PetscOptionsHasName(PETSC_NULL_CHARACTER,&
'-grouped_varr',this%gvarr,ierr)
call PetscOptionsGetString(PETSC_NULL_CHARACTER,&
'-matrix_schedule',this%sched,flg,ierr)
call PetscOptionsGetString(PETSC_NULL_CHARACTER,&
'-output_format',this%format,flg,ierr)
call PetscOptionsGetString(PETSC_NULL_CHARACTER,&
'-project_id',this%projid,flg,ierr)
call PetscOptionsGetString(PETSC_NULL_CHARACTER,&
'-inparam_file',this%inpar,flg,ierr)
call PetscOptionsGetString(PETSC_NULL_CHARACTER,&
'-synthetic_model',this%synth,flg,ierr)
call PetscOptionsGetString(PETSC_NULL_CHARACTER,&
'-varr_model_array',this%modarr,flg,ierr)
call PetscOptionsGetString(PETSC_NULL_CHARACTER,&
'-grid_info',this%grdin,flg,ierr)
call PetscOptionsGetString(PETSC_NULL_CHARACTER,&
'-solution_type',this%type,flg,ierr)
end subroutine parse_options
!========================================================
!
! identifies with which parameterization we are dealing with
!
subroutine tokenize_param(this)
implicit none
class(opts) :: this
PetscInt pos1
PetscInt pos2
PetscInt i
pos1 = 1
this%npars = 0
do
pos2 = index(this%param(pos1:), ",")
if (pos2 == 0) then
this%npars = this%npars + 1
this%parameter(this%npars) = this%param(pos1:)
exit
end if
this%npars = this%npars + 1
this%parameter(this%npars) = this%param(pos1:pos1+pos2-2)
pos1 = pos2+pos1
end do
end subroutine tokenize_param
!========================================================
end module module_options
!============================================================
!============================================================
!
! The module containing everything related to the mesh
!
module module_mesh
use global_param
use module_options
implicit none
type,public :: mesh
PetscInt nlays
PetscScalar eqinc
PetscBool adapt
PetscChar(256) grdin
PetscScalar eqinc_ref
PetscInt nlatzones
! To be defined
PetscInt blocks_per_param ! number of blocks per parameter
PetscInt blocks_all_param ! number of blocks per parameter
PetscInt, allocatable :: blocks_per_layer(:) ! number of blocks in each layer
PetscInt, allocatable :: nsqrs(:)
PetscInt, allocatable :: nsqrs_tot(:)
PetscScalar, allocatable :: xlamin(:,:)
PetscScalar, allocatable :: xlamax(:,:)
PetscScalar, allocatable :: xlomin(:,:)
PetscScalar, allocatable :: xlomax(:,:)
PetscScalar, allocatable :: locent(:,:)
PetscScalar, allocatable :: lacent(:,:)
PetscScalar, allocatable :: radmin(:,:)
PetscScalar, allocatable :: radmax(:,:)
PetscInt, allocatable :: levels(:,:)
contains
procedure :: setup_mesh
procedure :: read_mesh
procedure :: gen_mesh
procedure :: coordinates => get_coordinates
end type mesh
contains
!========================================================
!
! procedure to initialize a new mesh
!
subroutine setup_mesh(this,inopts)
implicit none
class(mesh) :: this
class(opts) :: inopts
! Local variables
PetscBool exist
PetscInt i,j
! Allocate memory
allocate ( this%blocks_per_layer( inopts%nlays ) )
allocate ( this%levels( n0max , inopts%nlays ) )
allocate ( this%xlamin( n0max , inopts%nlays ) )
allocate ( this%xlamax( n0max , inopts%nlays ) )
allocate ( this%xlomin( n0max , inopts%nlays ) )
allocate ( this%xlomax( n0max , inopts%nlays ) )
allocate ( this%radmin( n0max , inopts%nlays ) )
allocate ( this%radmax( n0max , inopts%nlays ) )
allocate ( this%locent( n0max , inopts%nlays ) )
allocate ( this%lacent( n0max , inopts%nlays ) )
allocate ( this%nsqrs_tot( nlmax ) ) ! 180/0.625
allocate ( this%nsqrs( nlmax ) ) ! 180/0.625
! Initialize more mesh params
this%eqinc = inopts%eqinc
this%adapt = inopts%adapt
this%nlays = inopts%nlays
this%eqinc_ref = inopts%eqinc_ref
this%blocks_per_layer(:) = 0
this%blocks_per_param = 0
! Select between variable and regular mesh
select case ( inopts%adapt )
case (.true.)
call PetscPrintf(PETSC_COMM_WORLD,&
"MESHER: Reading mesh from disk\n",ierr)
inquire ( file=inopts%grdin, exist=exist)
if (.not. exist) then
call PetscPrintf(PETSC_COMM_WORLD,&
"**** MESHER CRASHED **** : You didn't\
provide me with a grid info file!\n",ierr)
stop
end if
call this%read_mesh(inopts%grdin)
case (.false.)
! Generate a new grid based on eqincr
call PetscPrintf(PETSC_COMM_WORLD,&
"MESHER: Generating mesh on the fly!\n",ierr)
if ( inopts%eqinc.le.0.d5 ) then
call PetscPrintf(PETSC_COMM_WORLD,&
"**** MESHER CRASHED ****: You didn't\
provide me with a valid eqincr value!\n",ierr)
stop
end if
call this%gen_mesh()
end select
call PetscPrintf(PETSC_COMM_WORLD,"MESHER: Matrices are computed for "//&
trim(int2str(inopts%npars))//" parameters\n",ierr)
this%blocks_all_param = this%blocks_per_param*&
inopts%npars
end subroutine setup_mesh
!========================================================
!========================================================
!
! subroutine to read a mesh from ascii on disk
!
subroutine read_mesh(this, grdin)
implicit none
class(mesh) :: this
PetscChar(256) grdin
PetscChar(5) dummy
PetscInt i,j,ios
PetscInt lev1,lev2
PetscInt, parameter :: fz1=77,fz2=88,fz3=99 ! file handler
ios = 0
open(unit=fz1,file=trim(grdin),iostat=ios)
open(unit=fz2,file=trim(grdin)//".sh",iostat=ios)
open(unit=fz3,file=trim(grdin)//".lay",iostat=ios)
! Read a stored grid from disk
do i = 1,this%nlays
read(fz3,fmt=*), this%blocks_per_layer(i)
this%blocks_per_param = this%blocks_per_param + &
this%blocks_per_layer(i)
do j = 1,this%blocks_per_layer(i)
this%levels(j,i)=1
read(fz1, fmt=*, iostat=ios) this%levels(j,i)
read(fz2,fmt=*),dummy,this%xlamin(j,i),&
this%xlamax(j,i),&
this%xlomin(j,i),&
this%xlomax(j,i)
! Compute center lat and lon
this%lacent(j,i)=abs((this%xlamax(j,i)-&
this%xlamin(j,i))/2) + this%xlamin(j,i)
this%locent(j,i)=abs((this%xlomax(j,i)-&
this%xlomin(j,i))/2) + this%xlomin(j,i)
end do
end do
close(fz1)
close(fz2)
close(fz3)
call PetscPrintf(PETSC_COMM_WORLD,&
"MESHER: Mesh read successfully!\n",ierr)
end subroutine read_mesh
!========================================================
!========================================================
!
! subroutine to generate mesh on the fly
!
subroutine gen_mesh(this) !eqinc, eqinc_ref, nlays)
implicit none
class(mesh) :: this
PetscScalar coords( 6 )
PetscInt nsqrs_ref( nlmax )
PetscScalar colat, theta
PetscScalar delon, fact
PetscInt nlatz, nlatz_ref
PetscInt k, k_ref, n1lay
PetscInt i, j
nlatz_ref = int( 180. / this%eqinc_ref )
nlatz = int( 180. / this%eqinc )
colat = -this%eqinc_ref / 2.
! loop over latitudinal zones
do k = 1,nlatz_ref
colat = colat + this%eqinc_ref
theta = ( colat/180. )*pi
delon = this%eqinc_ref / ( sin(theta) )
nsqrs_ref(k) = int( 360./delon ) + 1
if (mod(nsqrs_ref(k),2).ne.0) nsqrs_ref(k) = nsqrs_ref(k)-1
! correct nsqrs(k) such to be compatible with reference grid
if (360./nsqrs_ref(k).ge.this%eqinc_ref) then
do while (mod(360./nsqrs_ref(k),this%eqinc_ref).ne.0)
nsqrs_ref(k) = nsqrs_ref(k) + 1
end do
elseif (360./nsqrs_ref(k).lt.this%eqinc_ref) then
do while (mod(this%eqinc_ref,360./nsqrs_ref(k)).ne.0)
nsqrs_ref(k) = nsqrs_ref(k) - 1
end do
endif
if (mod(nsqrs_ref(k),2).ne.0) then
call PetscPrintf(PETSC_COMM_WORLD,&
"**** MESHER CRASHED ****: # pixels per\
latitude has to be even\n",ierr)
stop
end if
enddo
! define actual grid as integer of reference grid
fact = this%eqinc_ref / this%eqinc
n1lay = 0
call PetscPrintf(PETSC_COMM_WORLD,&
"MESHER: Ratio of coarse to fine grid: "//&
trim(int2str(int(fact)))//"\n",ierr)
do k = 1,nlatz
k_ref = ( (k-1) / int(fact) ) + 1
this%nsqrs(k) = nsqrs_ref(k_ref) * int(fact)
n1lay = n1lay + this%nsqrs(k)
this%nsqrs_tot(k) = n1lay
enddo
this%nlatzones = nlatz
! in all layers we have the same nr of blocks
this%blocks_per_layer(:) = n1lay
call PetscPrintf(PETSC_COMM_WORLD,&
'MESHER: Number of pixel with finest\
parameterization: '//trim(int2str(&
this%blocks_per_layer(1)))//'\n',ierr)
! continue down to other layers
do i=1,this%nlays
do j=1,this%blocks_per_layer(i)
! Get coordinate ranges
coords = this%coordinates(j)
! Minimum and maximum lat and lon
this%levels(j,i) = 1 ! we just have one level
this%xlamin(j,i) = coords(1)
this%xlamax(j,i) = coords(2)
this%xlomin(j,i) = coords(3)
this%xlomax(j,i) = coords(4) ! dont need coords(5:6)
! Compute center lat and lon
this%lacent(j,i) = coords(5)
this%locent(j,i) = coords(6)
end do
this%blocks_per_param = this%blocks_per_param + &
this%blocks_per_layer(i)
end do
end subroutine gen_mesh
!========================================================
!========================================================
!
! procedure to get center and range of a block
!
function get_coordinates(this,isqr) result(coords)
implicit none
class(mesh), intent(in) :: this
PetscInt, intent(in) :: isqr
PetscScalar coords(6)
PetscScalar xlamax,xlamin
PetscScalar xlomax,xlomin
PetscScalar xlomid,xlamid
PetscScalar rlong,rlati,rinlo
PetscInt ila,isq,ntot
ntot = 0
xlomid = 0
xlomax = 0
xlomin = 0
xlamid = 0
xlamax = 0
xlamin = 0
! loop(s) over all the blocks
do ila=1,this%nlatzones
! increment latitude
rlati=90.-(this%eqinc*(ila-1))
! calculate increment in longitude for this band
rinlo=(360./this%nsqrs(ila))
do isq=1,this%nsqrs(ila)
rlong=(360./this%nsqrs(ila))*(isq-1)
ntot=ntot+1
if (ntot.eq.isqr) then
xlomid=rlong+(rinlo/2.d0)
xlamid=rlati-(this%eqinc/2.d0)
xlomax=rlong+(rinlo)
xlamax=rlati
xlomin=rlong
xlamin=rlati-(this%eqinc)
goto 600
end if
end do
end do
600 continue ! gets here in case of top or bottom layer
coords(1) = xlamin ! min lat
coords(2) = xlamax ! max lat
coords(3) = xlomin ! min lon
coords(4) = xlomax ! max lon
coords(5) = xlamid ! center lat
coords(6) = xlomid ! center lon
end function get_coordinates
!========================================================
end module module_mesh
!============================================================
!============================================================
!
! module related to the inversion schedule
!
module module_schedule
use module_mesh
use module_options
use global_param
implicit none
type,public :: sche
PetscInt mats_total ! number of submatrices
PetscInt rows_total ! total number of rows
PetscInt loop_total
PetscInt nrdamprows
PetscBool irdamp
PetscBool indamp
PetscBool iddamp
PetscBool iscale
PetscBool i3dref
PetscBool isynth
PetscChar(256) major_mode
PetscChar(256) inpar_file
PetscChar(256) sched_file
PetscScalar, allocatable :: ndamp(:,:)
PetscScalar, allocatable :: rdamp(:,:)
PetscScalar, allocatable :: ddamp(:,:)
PetscScalar, allocatable :: scale(:,:)
PetscScalar, allocatable :: layer(:)
PetscScalar, allocatable :: dloop(:,:)
PetscInt, allocatable :: prealloc(:)
PetscScalar, allocatable :: meta_info(:,:)
character(500), allocatable :: path_info(:,:)
PetscInt, allocatable :: fromto_info(:,:)
contains
procedure :: initialize_inversion
procedure :: read_schedule_file
procedure :: read_inparam_file
procedure :: memory_prealloc
end type sche
contains
!========================================================
!
! procededure to initialize all inversion schedules
!
subroutine initialize_inversion(this,inopts,inmesh)
implicit none
class(sche) :: this
class(opts) :: inopts
class(mesh) :: inmesh
this%inpar_file = inopts%inpar
this%sched_file = inopts%sched
this%loop_total = 0 ! number of global damping weights
this%rows_total = 0 ! total number of measurements
this%mats_total = 0 ! total number of submatrices
this%nrdamprows = 0 ! default is 0
! logicals on the type of damping applied
this%irdamp = .false.
this%indamp = .false.
this%iddamp = .false.
this%iscale = .false.
this%i3dref = .false.
! select major mode
this%major_mode = 'INVERSION'
if (trim(inopts%modarr).ne.'') then
this%major_mode = 'MODELARRAY'
end if
! logicals about synthetic testing
if (trim(inopts%synth)=='') then
this%isynth = .false.
else
this%isynth = .true.
end if
! is the reference model 3D?
inquire(file=trim(inopts%refmod)//"/regions.dat", exist=this%i3dref)
! allocate submatrix and damping schedule
allocate ( this%rdamp( inopts%nlays, npmax+1 ) )
allocate ( this%ndamp( inopts%nlays, npmax ) )
allocate ( this%ddamp( inopts%nlays, 2 ) )
allocate ( this%scale( inopts%nlays, 3 ) )
allocate ( this%layer( inopts%nlays+1 ) )
allocate ( this%path_info( matmx, 4 ) )
allocate ( this%meta_info( matmx, 4 ) )
allocate ( this%fromto_info( matmx, 2 ) )
allocate ( this%dloop( 4, dmpmx ) )
! read damping and matrix schedules
call this % read_inparam_file ( inopts, inmesh )
call this % read_schedule_file ( inopts )
! now we now how large the matrix is going to be
! so allocate space for preallocvec
allocate ( this%prealloc ( this%nrdamprows + &
this%rows_total ) )
call this % memory_prealloc ( )
end subroutine initialize_inversion
!========================================================
!========================================================
!
! Procededure to read a damping parameters and other
! necessary information from the inparam file
!
subroutine read_inparam_file(this,inopts,inmesh)
implicit none
class(sche):: this
class(opts):: inopts
class(mesh):: inmesh
PetscInt ios,j,i
PetscInt, parameter :: fh=16 ! file handler
PetscChar(256) line,key,val,dummy
open(unit=fh,file=this%inpar_file,iostat=ios)
ios = 0
do
read(fh, fmt='(a256)', iostat=ios) line
if (ios < 0) then
call PetscPrintf(PETSC_COMM_WORLD,&
"SCHEDULER: Read damping scheme successfully\n",ierr)
exit
end if
if (len(trim(line)) < 1 .or.&
line(1:1) == ';') then
continue
else
select case(trim(line))
case ('RDAMP_DEPTH')
call PetscPrintf(PETSC_COMM_WORLD,&
"SCHEDULER: Reading roughness damping scheme\n",ierr)
read(fh, fmt=*, iostat=ios) dummy
do i=1,inopts%nlays
read(fh, fmt=*, iostat=ios) this%rdamp(i,1:5)
if ( verbosity > 1 ) print*,this%rdamp(i,:)
end do
case ('NDAMP_DEPTH')
call PetscPrintf(PETSC_COMM_WORLD,&
"SCHEDULER: Reading norm damping scheme\n",ierr)
read(fh, fmt=*, iostat=ios) dummy
do i=1,inopts%nlays
read(fh, fmt=*, iostat=ios) this%ndamp(i,1:4)
if ( verbosity > 1 ) print*,this%ndamp(i,:)
end do
case ('DDAMP_DEPTH')
call PetscPrintf(PETSC_COMM_WORLD,&
"SCHEDULER: Reading difference damping scheme\n",ierr)
read(fh, fmt=*, iostat=ios) dummy
do i=1,inopts%nlays
read(fh, fmt=*, iostat=ios) this%ddamp(i,1:2)
if ( verbosity > 1 ) print*,this%ddamp(i,1:2)
end do
case ('SCALING_COEFFS')
call PetscPrintf(PETSC_COMM_WORLD,&
"SCHEDULER: Reading vp-to-vs scaling scheme\n",ierr)
read(fh, fmt=*, iostat=ios) dummy
do i=1,inopts%nlays
read(fh, fmt=*, iostat=ios) this%scale(i,1:3)
if ( verbosity > 1 ) print*,this%scale(i,1:3)
end do
case ('LAYER_BOTTOMS')
call PetscPrintf(PETSC_COMM_WORLD,&
"SCHEDULER: Reading layer bottom depths for postprocessing\n",ierr)
read(fh, fmt=*, iostat=ios) dummy
this%layer(1)=0.d0 ! 0th index interface always surface = depth zero
do i=2,inopts%nlays+1
read(fh, fmt=*, iostat=ios) this%layer(i)
if ( verbosity > 1 ) print*,this%layer(i)
end do
case ('DAMP_LOOP')
call PetscPrintf(PETSC_COMM_WORLD,&
"SCHEDULER: Reading global damping loop scheme\n",ierr)
read(fh, fmt=*, iostat=ios) dummy
read(fh, fmt=*, iostat=ios) this%loop_total
do i=1,4
read(fh, fmt=*, iostat=ios) &
this%dloop(i,1:this%loop_total)
if ( verbosity > 1 ) print*,this%dloop(i,1:this%loop_total)
end do
! determine total number of damping rows
do i=1,this%loop_total
! nr of roughness damping rows
if ( this%dloop(1,i) > 0. ) this%irdamp = .true.
if ( this%dloop(2,i) > 0. ) this%indamp = .true.
if ( this%dloop(3,i) > 0. ) this%iddamp = .true.
if ( this%dloop(4,i) > 0. ) this%iscale = .true.
end do
if ( this%indamp ) then
this%nrdamprows = this%nrdamprows + &
inmesh%blocks_per_param * &
inopts%npars
end if
if ( this%irdamp ) then
this%nrdamprows = this%nrdamprows + &
inmesh%blocks_per_param * &
inopts%npars
end if
if ( this%iddamp ) then
if (inopts%npars==4) then
this%nrdamprows = this%nrdamprows + &
inmesh%blocks_per_param*2
else if ((inopts%npars==2).or.(inopts%npars==3)) then
this%nrdamprows = this%nrdamprows + &
inmesh%blocks_per_param
end if
end if
if ( this%iscale ) then
if (inopts%npars==4) then
this%nrdamprows = this%nrdamprows + &
inmesh%blocks_per_param*2
else if (inopts%npars==2) then
this%nrdamprows = this%nrdamprows + &
inmesh%blocks_per_param
end if
end if
call PetscPrintf(PETSC_COMM_WORLD,"SCHEDULER: In total we attach : "//&
trim(int2str(this%nrdamprows))//" damping rows\n",ierr)
end select
end if
end do
close(fh)
end subroutine read_inparam_file
!========================================================
!========================================================
!
! procededure to read a matrix schedule
!
subroutine read_schedule_file(this,inopts)
implicit none
class(sche):: this
class(opts):: inopts
PetscInt ios,j,i
PetscInt, parameter :: fh=15 ! file handler
ios = 0
open(fh,file=trim(this%sched_file))
this%mats_total = 1
do
call PetscPrintf(PETSC_COMM_WORLD,&
"SCHEDULER: Reading meta info for submatrix "//&
trim(int2str(this%mats_total))//"\n",ierr)
do j=1,4
read(unit=fh,fmt=*,iostat=ios) &
this%path_info( this%mats_total, j )
if (verbosity >1 ) print*,trim(this%path_info( this%mats_total, j ))
end do
do j=1,4
read(unit=fh,fmt=*,iostat=ios) &
this%meta_info( this%mats_total, j )
if (verbosity > 1 ) print*,this%meta_info( this%mats_total, j )
end do
if ( trim(this%path_info( this%mats_total, 1 )) =='EOF' ) then
! remove trash written in the last entry
this%meta_info( this%mats_total, : ) = 0.
this%path_info( this%mats_total, : ) = ''
this%mats_total = this%mats_total - 1
exit
else
call PetscPrintf(PETSC_COMM_WORLD,"SCHEDULER: Subset contains "//&
trim(int2str(int(this%meta_info(this%mats_total,4))))//" rows\n",ierr)
this%rows_total = this%rows_total + &
int(this%meta_info( this%mats_total, 4 ))
this%mats_total = this%mats_total + 1
end if
end do
close(fh)
call PetscPrintf(PETSC_COMM_WORLD,&
"SCHEDULER: We will combine a total of: "//&
trim(int2str(this%rows_total))//" data\n",ierr)
call PetscPrintf(PETSC_COMM_WORLD,&
"SCHEDULER: Successfully read submatrix schedule!\n ",ierr)
end subroutine read_schedule_file
!========================================================
!========================================================
!
! procededure to define
!
subroutine memory_prealloc(this)
implicit none
class(sche):: this
PetscInt i,j,row_glob
PetscInt, allocatable :: pnt(:)
! loop over all submatrices, detailed preallocation
row_glob = 0
do i=1,this%mats_total
! meta_info(i,3) contains nr rows of submatrix
allocate( pnt ( int(this%meta_info(i,4)) ) )
! path_info(i,3) contains pointer vector
open(70,file=trim(this%path_info(i,3)),status='old')
! loop over rows in submatrix
do j=1,int(this%meta_info(i,4))
row_glob=row_glob+1
read(70,*) pnt(j)
this%prealloc(row_glob) = pnt(j)-pnt(j-1)
if (this%prealloc(row_glob).gt.nprmx) then
print*,"PREALLOC CRASHED: nprmx too small"
stop
end if
end do
deallocate ( pnt )
close(70)
end do
! and now the same for all damping rows
do i=1,this%nrdamprows
row_glob=row_glob+1
! max nr of entries per damping row
this%prealloc(row_glob) = nperd
end do
call PetscPrintf(PETSC_COMM_WORLD,"SCHEDULER: Successfully"//&
"pre-allocated memory!\n ",ierr)
end subroutine memory_prealloc
!========================================================
end module module_schedule
!============================================================
!============================================================
!
! module containing everything related to reading matrices
!
module module_matrix
use global_param
use module_schedule
use module_options
use module_mesh
implicit none
type,public :: matr
Mat A ! global system matrix
Vec x ! solution vector
Vec b ! rhs vector
Vec b_store ! to store the rhs vector
Vec b_synth ! synthetics rhs vector
Vec b_adotx ! rhs vector for cummulative variance red
Vec b_dummy ! rhs vector for grouped variance red
Vec b_rough ! synthetics rhs vector
Vec x_synth ! synthetics model vector
Vec mask_dmp ! masking the damping rows
Vec mask_sub ! masking everything but subset rows
Vec mask_dat ! masking everything but roughness damping rows
Mat ATb ! A transposed * b
Mat ATA ! bormal equations
PetscInt row
PetscInt row_start
PetscInt row_end
PetscMPIInt processor
contains
procedure :: initialize_matrix
procedure :: read_submatrices
procedure :: assemble_matrix
procedure :: assemble_vectors
procedure :: apply_rdamp
procedure :: apply_ddamp
procedure :: apply_ndamp
procedure :: apply_scale
procedure :: store_rhs
procedure :: restore_rhs
procedure :: compute_adotx
procedure :: destroy_matrix
procedure :: read_synth_model
end type matr
contains
!========================================================
!
! procededure to initialize linear system to solve
!
subroutine initialize_matrix(this,insche,inmesh)
implicit none
class(matr) :: this
class(sche) :: insche
class(mesh) :: inmesh
! Get the rank of this processor
call MPI_Comm_rank ( PETSC_COMM_WORLD, this%processor, ierr)
call MatCreate ( PETSC_COMM_WORLD, this%A, ierr )
! changing size of A matrix, to make space for damping rows
call MatSetSizes ( this%A, PETSC_DECIDE, &
PETSC_DECIDE, insche%rows_total+ &
insche%nrdamprows, &
inmesh%blocks_all_param, ierr)