-
Notifications
You must be signed in to change notification settings - Fork 119
Expand file tree
/
Copy pathSolvers.F90
More file actions
2874 lines (2412 loc) · 105 KB
/
Solvers.F90
File metadata and controls
2874 lines (2412 loc) · 105 KB
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 (C) 2006 Imperial College London and others.
!
! Please see the AUTHORS file in the main source directory for a full list
! of copyright holders.
!
! Prof. C Pain
! Applied Modelling and Computation Group
! Department of Earth Science and Engineering
! Imperial College London
!
! amcgsoftware@imperial.ac.uk
!
! This library is free software; you can redistribute it and/or
! modify it under the terms of the GNU Lesser General Public
! License as published by the Free Software Foundation; either
! version 2.1 of the License, or (at your option) any later version.
!
! This library is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
! Lesser General Public License for more details.
!
! You should have received a copy of the GNU Lesser General Public
! License along with this library; if not, write to the Free Software
! Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
! USA
#include "fdebug.h"
module solvers
use FLDebug
use Global_Parameters
use futils, only: present_and_true, int2str, free_unit, real_format
use element_numbering, only: ELEMENT_BUBBLE
use elements
use spud
use parallel_tools
use petsc
use Sparse_Tools
use fields_calculations
use Fields
use profiler
use Petsc_tools
use Signal_Vars
use Multigrid
use sparse_tools_petsc
use sparse_matrices_fields
use vtk_interfaces
use halos
use MeshDiagnostics
implicit none
! Module to provide explicit interfaces to matrix solvers.
#include "petsc_legacy.h"
! stuff used in the PETSc monitor (see petsc_solve_callback_setup() below)
integer :: petsc_monitor_iteration = 0
Vec :: petsc_monitor_x
!
! if .true. the code will compare with the provided exact answer, and
! give the error convergence each iteration:
logical, save:: petsc_monitor_has_exact=.false.
! this requires the following:
Vec :: petsc_monitor_exact
real, dimension(:), pointer :: petsc_monitor_error => null()
PetscLogDouble, dimension(:), pointer :: petsc_monitor_flops => null()
type(scalar_field), save:: petsc_monitor_exact_sfield
type(vector_field), save:: petsc_monitor_exact_vfield
character(len=FIELD_NAME_LEN), save:: petsc_monitor_error_filename=""
!
! if .true. a vtu will be written for each iteration
logical, save:: petsc_monitor_iteration_vtus=.false.
! this requires the following:
type(petsc_numbering_type), save:: petsc_monitor_numbering
type(vector_field), target, save:: petsc_monitor_positions
type(scalar_field), dimension(3), save:: petsc_monitor_sfields
type(vector_field), dimension(3), save:: petsc_monitor_vfields
character(len=FIELD_NAME_LEN), save:: petsc_monitor_vtu_name
integer, save:: petsc_monitor_vtu_series=0
private
public petsc_solve, set_solver_options, &
complete_solver_option_path, petsc_solve_needs_positions, &
L2_project_nullspace_vector
! meant for unit-testing solver code only:
public petsc_solve_core, petsc_solve_destroy, &
petsc_solve_copy_vectors_from_scalar_fields, &
setup_ksp_from_options, create_ksp_from_options, petsc_solve_monitor_exact, &
petsc_solve_monitor_iteration_vtus, attach_null_space_from_options, &
petsc_solve_setup
interface petsc_solve
module procedure petsc_solve_scalar, petsc_solve_vector, &
petsc_solve_scalar_multiple, &
petsc_solve_vector_components, &
petsc_solve_tensor_components, &
petsc_solve_scalar_petsc_csr, petsc_solve_vector_petsc_csr
end interface
interface set_solver_options
module procedure set_solver_options_with_path, &
set_solver_options_scalar, set_solver_options_vector, set_solver_options_tensor
end interface set_solver_options
interface petsc_solve_monitor_exact
module procedure petsc_solve_monitor_exact_scalar
end interface petsc_solve_monitor_exact
contains
subroutine petsc_solve_scalar(x, matrix, rhs, option_path, &
preconditioner_matrix, prolongators, surface_node_list, &
internal_smoothing_option, iterations_taken)
!!< Solve a linear system the nice way.
type(scalar_field), intent(inout) :: x
type(scalar_field), intent(in) :: rhs
type(csr_matrix), intent(in) :: matrix
character(len=*), optional, intent(in) :: option_path
!! 2 experimental arguments to improve preconditioning with extra outside information
!! provide approximation the matrix (only to be used in combination with pctype='KSP')
type(csr_matrix), optional, intent(in) :: preconditioner_matrix
!! prolongators to be used at the first levels of 'mg'
type(petsc_csr_matrix), dimension(:), optional, intent(in) :: prolongators
!! surface_node_list for internal smoothing
integer, dimension(:), optional, intent(in) :: surface_node_list
!! internal smoothing option
integer, intent(in), optional :: internal_smoothing_option
!! the number of petsc iterations taken
integer, intent(out), optional :: iterations_taken
KSP ksp
Mat A
Vec y, b
character(len=OPTION_PATH_LEN):: solver_option_path
type(petsc_numbering_type) petsc_numbering
integer literations
logical lstartfromzero
assert(size(x%val)==size(rhs%val))
assert(size(x%val)==size(matrix,2))
#ifdef DDEBUG
if (.not.associated(matrix%sparsity%column_halo)) then
assert(size(rhs%val)==size(matrix,1))
else
! in parallel we allow for the matrix to not have the halo nodes
! whereas the rhs will always contain them (even if not used)
if (size(matrix,1)/=size(rhs%val)) then
! get_nowned_nodes() seems completely buggered for T10 halo
! this should be fixed by dham+jrmaddison's mesh+halo integration
!assert(size(matrix,1)==get_nowned_nodes(matrix%sparsity%halo_tag))
assert( size(matrix,1)==halo_nowned_nodes(matrix%sparsity%column_halo) )
end if
end if
#endif
! setup PETSc object and petsc_numbering from options and
call petsc_solve_setup(y, A, b, ksp, petsc_numbering, &
solver_option_path, lstartfromzero, &
matrix=matrix, &
sfield=x, &
option_path=option_path, &
preconditioner_matrix=preconditioner_matrix, &
prolongators=prolongators, surface_node_list=surface_node_list, &
internal_smoothing_option=internal_smoothing_option)
! copy array into PETSc vecs
call petsc_solve_copy_vectors_from_scalar_fields(y, b, x, &
& matrix, rhs, petsc_numbering, lstartfromzero)
! the solve and convergence check
call petsc_solve_core(y, A, b, ksp, petsc_numbering, &
solver_option_path, lstartfromzero, literations, &
sfield=x, x0=x%val)
! set the optional variable passed out of this procedure
! for the number of petsc iterations taken
if (present(iterations_taken)) iterations_taken = literations
! Copy back the result using the petsc numbering:
call petsc2field(y, petsc_numbering, x, rhs)
! destroy all PETSc objects and the petsc_numbering
call petsc_solve_destroy(y, A, b, ksp, petsc_numbering, &
& solver_option_path)
end subroutine petsc_solve_scalar
subroutine petsc_solve_scalar_multiple(x, matrix, rhs, option_path)
!!< Solves multiple scalar fields with the same matrix.
!!< Need to specify an option_path as there's no default
type(scalar_field), dimension(:), intent(inout) :: x
type(scalar_field), dimension(:), intent(in) :: rhs
type(csr_matrix), intent(in) :: matrix
character(len=*), optional, intent(in) :: option_path
KSP ksp
Mat A
Vec y, b
type(petsc_numbering_type) petsc_numbering
character(len=OPTION_PATH_LEN) solver_option_path
integer literations
logical lstartfromzero
integer i
assert(size(x)==size(rhs))
do i=1, size(x)
assert(size(x(i)%val)==size(rhs(i)%val))
assert(size(x(i)%val)==size(matrix,2))
assert(size(rhs(i)%val)==size(matrix,1))
end do
ewrite(1,*) 'Solving for multiple scalar fields at once'
! setup PETSc object and petsc_numbering from options and
call petsc_solve_setup(y, A, b, ksp, petsc_numbering, &
solver_option_path, lstartfromzero, &
matrix=matrix, sfield=x(1), &
option_path=option_path)
do i=1, size(x)
! copy array into PETSc vecs
call petsc_solve_copy_vectors_from_scalar_fields(y, b, &
x(i), matrix, rhs(i), &
petsc_numbering, lstartfromzero)
! the solve and convergence check
call petsc_solve_core(y, A, b, ksp, petsc_numbering, &
solver_option_path, lstartfromzero, literations, &
sfield=x(i), x0=x(i)%val)
! Copy back the result using the petsc numbering:
call petsc2field(y, petsc_numbering, x(i), rhs(i))
end do
ewrite(1,*) 'Finished solving all scalar fields'
! destroy all PETSc objects and the petsc_numbering
call petsc_solve_destroy(y, A, b, ksp, petsc_numbering, &
& solver_option_path)
end subroutine petsc_solve_scalar_multiple
subroutine petsc_solve_vector(x, matrix, rhs, option_path, deallocate_matrix)
!!< Solve a linear system the nice way. Options for this
!!< come via the options mechanism.
type(vector_field), intent(inout) :: x
type(vector_field), intent(in) :: rhs
type(block_csr_matrix), intent(inout) :: matrix
character(len=*), optional, intent(in) :: option_path
!! deallocate the matrix after it's been copied
logical, intent(in), optional :: deallocate_matrix
KSP ksp
Mat A
Vec y, b
type(petsc_numbering_type) petsc_numbering
character(len=OPTION_PATH_LEN) solver_option_path
integer literations
logical lstartfromzero
type(csr_matrix) :: matrixblock
type(scalar_field) :: rhsblock, xblock
integer :: i
assert(x%dim==rhs%dim)
assert(size(x%val(1,:))==size(rhs%val(1,:)))
assert(size(x%val(1,:))==block_size(matrix,2))
assert(size(rhs%val(1,:))==block_size(matrix,1))
assert(x%dim==blocks(matrix,2))
assert(rhs%dim==blocks(matrix,1))
if(matrix%diagonal) then
assert(blocks(matrix,1)==blocks(matrix,2))
! only want to solve using the diagonal blocks
do i = 1, blocks(matrix,1)
matrixblock=block(matrix,i,i)
rhsblock = extract_scalar_field(rhs, i)
xblock = extract_scalar_field(x, i)
! setup PETSc object and petsc_numbering from options and
call petsc_solve_setup(y, A, b, ksp, petsc_numbering, &
solver_option_path, lstartfromzero, &
matrix=matrixblock, &
vfield=x, &
option_path=option_path)
! copy array into PETSc vecs
call petsc_solve_copy_vectors_from_scalar_fields(y, b, xblock, &
& matrixblock, rhsblock, petsc_numbering, lstartfromzero)
if(present_and_true(deallocate_matrix).and.(i==blocks(matrix,1))) then
call deallocate(matrix)
end if
! the solve and convergence check
call petsc_solve_core(y, A, b, ksp, petsc_numbering, &
solver_option_path, lstartfromzero, literations, &
vfield=x, x0=xblock%val)
! Copy back the result using the petsc numbering:
call petsc2field(y, petsc_numbering, xblock, rhsblock)
! destroy all PETSc objects and the petsc_numbering
call petsc_solve_destroy(y, A, b, ksp, petsc_numbering, &
& solver_option_path)
end do
else
! setup PETSc object and petsc_numbering from options and
call petsc_solve_setup(y, A, b, ksp, petsc_numbering, &
solver_option_path, lstartfromzero, &
block_matrix=matrix, &
vfield=x, &
option_path=option_path)
if(present_and_true(deallocate_matrix)) then
call deallocate(matrix)
end if
! copy array into PETSc vecs
call petsc_solve_copy_vectors_from_vector_fields(y, b, x, rhs, petsc_numbering, lstartfromzero)
! the solve and convergence check
call petsc_solve_core(y, A, b, ksp, petsc_numbering, &
solver_option_path, lstartfromzero, literations, &
vfield=x, vector_x0=x)
! Copy back the result using the petsc numbering:
call petsc2field(y, petsc_numbering, x)
! destroy all PETSc objects and the petsc_numbering
call petsc_solve_destroy(y, A, b, ksp, petsc_numbering, &
& solver_option_path)
end if
end subroutine petsc_solve_vector
subroutine petsc_solve_vector_components(x, matrix, rhs, option_path)
!!< Solve a linear system the nice way. Options for this
!!< come via the options mechanism. This version solves a linear system
!!< for each of the components of rhs each time with the same matrix.
type(vector_field), intent(inout) :: x
type(vector_field), intent(in) :: rhs
type(csr_matrix), intent(in) :: matrix
character(len=*), optional, intent(in) :: option_path
KSP ksp
Mat A
Vec y, b
type(scalar_field) x_component, rhs_component
type(petsc_numbering_type) petsc_numbering
character(len=OPTION_PATH_LEN) solver_option_path, option_path_in
integer literations, i
logical lstartfromzero
assert(x%dim==rhs%dim)
assert(size(x%val(1,:))==size(rhs%val(1,:)))
assert(size(x%val(1,:))==size(matrix,2))
assert(size(rhs%val(1,:))==size(matrix,1))
! option_path_in may still point to field
! (so we have to add "/prognostic/solver" below)
if (present(option_path)) then
option_path_in=option_path
else
option_path_in=x%option_path
end if
! setup PETSc object and petsc_numbering from options and
call petsc_solve_setup(y, A, b, ksp, petsc_numbering, &
solver_option_path, lstartfromzero, &
matrix=matrix, &
vfield=x, &
option_path=option_path)
ewrite(1,*) 'Solving for multiple components of a vector field'
do i=1, x%dim
ewrite(1, *) 'Now solving for component: ', i
x_component=extract_scalar_field(x, i)
rhs_component=extract_scalar_field(rhs, i)
! copy array into PETSc vecs
call petsc_solve_copy_vectors_from_scalar_fields(y, b, x_component, matrix, rhs_component, petsc_numbering, lstartfromzero)
! the solve and convergence check
call petsc_solve_core(y, A, b, ksp, petsc_numbering, &
solver_option_path, lstartfromzero, literations, &
vfield=x, x0=x_component%val)
! Copy back the result using the petsc numbering:
call petsc2field(y, petsc_numbering, x_component, rhs_component)
end do
ewrite(1,*) 'Finished solving all components.'
! destroy all PETSc objects and the petsc_numbering
call petsc_solve_destroy(y, A, b, ksp, petsc_numbering, &
& solver_option_path)
end subroutine petsc_solve_vector_components
subroutine petsc_solve_scalar_petsc_csr(x, matrix, rhs, option_path, &
prolongators, surface_node_list)
!!< Solve a linear system the nice way. Options for this
!!< come via the options mechanism.
type(scalar_field), intent(inout) :: x
type(scalar_field), intent(in) :: rhs
type(petsc_csr_matrix), intent(inout) :: matrix
character(len=*), optional, intent(in) :: option_path
!! prolongators to be used at the first levels of 'mg'
type(petsc_csr_matrix), dimension(:), optional, intent(in) :: prolongators
!! surface_node_list for internal smoothing
integer, dimension(:), optional, intent(in) :: surface_node_list
Vec y, b
character(len=OPTION_PATH_LEN) solver_option_path
integer literations
logical lstartfromzero
assert(size(x%val)==size(rhs%val))
assert(size(x%val)==size(matrix,2))
assert(size(rhs%val)==size(matrix,1))
! setup PETSc object and petsc_numbering from options and
call petsc_solve_setup_petsc_csr(y, b, &
solver_option_path, lstartfromzero, &
matrix, &
sfield=x, &
option_path=option_path, &
prolongators=prolongators, surface_node_list=surface_node_list)
! copy array into PETSc vecs
call petsc_solve_copy_vectors_from_scalar_fields(y, b, x, rhs=rhs, &
petsc_numbering=matrix%row_numbering, startfromzero=lstartfromzero)
! the solve and convergence check
call petsc_solve_core(y, matrix%M, b, matrix%ksp, matrix%row_numbering, &
solver_option_path, lstartfromzero, literations, &
sfield=x, x0=x%val)
! Copy back the result using the petsc numbering:
call petsc2field(y, matrix%column_numbering, x)
! destroy all PETSc objects and the petsc_numbering
call petsc_solve_destroy_petsc_csr(y, b, solver_option_path)
end subroutine petsc_solve_scalar_petsc_csr
subroutine petsc_solve_vector_petsc_csr(x, matrix, rhs, option_path, &
prolongators, positions, rotation_matrix)
!!< Solve a linear system the nice way. Options for this
!!< come via the options mechanism.
type(vector_field), intent(inout) :: x
type(vector_field), intent(in) :: rhs
type(petsc_csr_matrix), intent(inout) :: matrix
character(len=*), optional, intent(in) :: option_path
!! prolongators to be used at the first levels of 'mg'
type(petsc_csr_matrix), dimension(:), optional, intent(in) :: prolongators
!! positions field is only used with remove_null_space/ or multigrid_near_null_space/ with rotational components
type(vector_field), intent(in), optional :: positions
!! with rotated bcs: matrix to transform from x,y,z aligned vectors to boundary aligned
Mat, intent(in), optional:: rotation_matrix
Vec y, b
character(len=OPTION_PATH_LEN) solver_option_path
integer literations
logical lstartfromzero
assert(x%dim==rhs%dim)
assert(size(x%val(1,:))==size(rhs%val(1,:)))
assert(size(x%val(1,:))==block_size(matrix,2))
assert(size(rhs%val(1,:))==block_size(matrix,1))
assert(x%dim==blocks(matrix,2))
assert(rhs%dim==blocks(matrix,1))
! setup PETSc object and petsc_numbering from options and
call petsc_solve_setup_petsc_csr(y, b, &
solver_option_path, lstartfromzero, &
matrix, vfield=x, option_path=option_path, &
prolongators=prolongators, &
positions=positions, rotation_matrix=rotation_matrix)
! copy array into PETSc vecs
call petsc_solve_copy_vectors_from_vector_fields(y, b, x, rhs, &
matrix%row_numbering, lstartfromzero)
! the solve and convergence check
call petsc_solve_core(y, matrix%M, b, matrix%ksp, matrix%row_numbering, &
solver_option_path, lstartfromzero, literations, &
vfield=x, vector_x0=x)
! Copy back the result using the petsc numbering:
call petsc2field(y, matrix%column_numbering, x)
! destroy all PETSc objects and the petsc_numbering
call petsc_solve_destroy_petsc_csr(y, b, solver_option_path)
end subroutine petsc_solve_vector_petsc_csr
subroutine petsc_solve_tensor_components(x, matrix, rhs, &
symmetric, option_path)
!!< Solve a linear system the nice way. Options for this
!!< come via the options mechanism. This version solves a linear system
!!< for each of the components of rhs each time with the same matrix.
type(tensor_field), intent(inout) :: x
type(tensor_field), intent(in) :: rhs
type(csr_matrix), intent(in) :: matrix
! if .true. assume rhs is symmetric (so we need to solve for fewer components)
logical, optional, intent(in):: symmetric
character(len=*), optional, intent(in) :: option_path
KSP ksp
Mat A
Vec y, b
type(scalar_field) x_component, rhs_component
type(petsc_numbering_type) petsc_numbering
character(len=OPTION_PATH_LEN) solver_option_path, option_path_in
integer literations, i, j, startj
logical lstartfromzero
assert(all(x%dim==rhs%dim))
assert(size(x%val,3)==size(rhs%val,3))
assert(size(x%val,3)==size(matrix,2))
assert(size(rhs%val,3)==size(matrix,1))
! option_path_in may still point to field
! (so we have to add "/prognostic/solver" below)
if (present(option_path)) then
option_path_in=option_path
else
option_path_in=x%option_path
end if
! setup PETSc object and petsc_numbering from options and
call petsc_solve_setup(y, A, b, ksp, petsc_numbering, &
solver_option_path, lstartfromzero, &
matrix=matrix, &
tfield=x, &
option_path=option_path_in)
ewrite(1,*) 'Solving for multiple components of a tensor field'
startj=1
do i=1, x%dim(1)
if (present(symmetric)) then
if (symmetric) then
! only computes with rhs(i,j) where j>=i
startj=i
end if
end if
do j=startj, x%dim(2)
ewrite(1, *) 'Now solving for component: ', i, j
x_component=extract_scalar_field(x, i, j)
rhs_component=extract_scalar_field(rhs, i, j)
! copy array into PETSc vecs
call petsc_solve_copy_vectors_from_scalar_fields(y, b, x_component, matrix, rhs_component, petsc_numbering, lstartfromzero)
! the solve and convergence check
call petsc_solve_core(y, A, b, ksp, petsc_numbering, &
solver_option_path, lstartfromzero, literations, &
tfield=x, x0=x_component%val)
! Copy back the result using the petsc numbering:
call petsc2field(y, petsc_numbering, x_component, rhs_component)
end do
end do
ewrite(1,*) 'Finished solving all components.'
if (present(symmetric)) then
if (symmetric) then
ewrite(2,*) 'This is a symmetric matrix'
ewrite(2,*) 'Only components (i,j) with j>=i have been solved for.'
ewrite(2,*) 'Now copying these to components (j,i).'
! copy results x(i,j) of equations with rhs(i,j) where j>=i to x(j,i)
do i=1, x%dim(1)
do j=i, x%dim(2)
x%val(j,i,:)=x%val(i,j,:)
end do
end do
end if
end if
! destroy all PETSc objects and the petsc_numbering
call petsc_solve_destroy(y, A, b, ksp, petsc_numbering, solver_option_path)
end subroutine petsc_solve_tensor_components
function complete_solver_option_path(option_path)
character(len=*), intent(in):: option_path
character(len=OPTION_PATH_LEN):: complete_solver_option_path
! at the moment only prognostic fields have a solver options block
! under [field_path]/prognostic/solver/. Other cases should be
! implemented here:
if (have_option(trim(option_path)//'/prognostic/solver')) then
complete_solver_option_path=trim(option_path)//'/prognostic/solver'
! some diagnostic cases have solver blocks now
else if (have_option(trim(option_path)//'/diagnostic/solver')) then
complete_solver_option_path=trim(option_path)//'/diagnostic/solver'
else if (have_option(trim(option_path)//'/solver')) then
complete_solver_option_path=trim(option_path)//'/solver'
else
ewrite(-1,*) 'option_path: ', trim(option_path)
FLAbort("Missing solver element in provided option_path.")
end if
end function complete_solver_option_path
subroutine petsc_solve_setup(y, A, b, ksp, petsc_numbering, &
solver_option_path, startfromzero, &
matrix, block_matrix, sfield, vfield, tfield, &
option_path, startfromzero_in, &
preconditioner_matrix, prolongators, surface_node_list, &
internal_smoothing_option, positions)
!!< sets up things needed to call petsc_solve_core
!! Stuff that comes out:
!!
!! PETSc solution vector
Vec, intent(out):: y
!! PETSc matrix
Mat, intent(out):: A
!! PETSc rhs vector
Vec, intent(out):: b
!! Solver object
KSP, intent(out):: ksp
!! numbering from local (i.e. fluidity speak: global) to PETSc (fluidity: universal) numbering
type(petsc_numbering_type), intent(out):: petsc_numbering
!! returns the option path to solver/ block for new options, otherwise ""
character(len=*), intent(out):: solver_option_path
!! whether to start with zero initial guess
logical, intent(out):: startfromzero
!! Stuff that goes in:
!!
!! provide either a matrix or block_matrix to be solved
type(csr_matrix), target, optional, intent(in):: matrix
type(block_csr_matrix), target, optional, intent(in):: block_matrix
!! provide either a scalar field or vector field to be solved for
type(scalar_field), optional, intent(in):: sfield
type(vector_field), optional, intent(in):: vfield
type(tensor_field), optional, intent(in):: tfield
!! if provided overrides sfield%option_path
character(len=*), optional, intent(in):: option_path
!! whether to start with zero initial guess (as passed in)
logical, optional, intent(in):: startfromzero_in
!! provide approximation the matrix (only to be used in combination with pctype='KSP')
type(csr_matrix), optional, intent(in) :: preconditioner_matrix
!! prolongators to be used at the first level of 'mg'
type(petsc_csr_matrix), dimension(:), optional, intent(in) :: prolongators
!! Stuff needed for internal smoother
integer, dimension(:), optional, intent(in) :: surface_node_list
integer, optional, intent(in) :: internal_smoothing_option
!! positions field is only used with remove_null_space/ or multigrid_near_null_space/ with rotational components
type(vector_field), intent(in), optional :: positions
logical, dimension(:), pointer:: inactive_mask
integer, dimension(:), allocatable:: ghost_nodes
Mat:: pmat
! one of the PETSc supplied orderings see
! http://www-unix.mcs.anl.gov/petsc/petsc-as/snapshots/petsc-current/docs/manualpages/MatOrderings/MatGetOrdering.html
MatOrderingType:: ordering_type
PetscBool:: use_reordering
real time1, time2
integer ierr
logical:: parallel, timing, have_cache
type(halo_type), pointer :: halo
integer i, j
character(len=FIELD_NAME_LEN) :: name
KSP, pointer:: ksp_pointer
! Initialise profiler
if(present(sfield)) then
call profiler_tic(sfield, "petsc_setup")
name = sfield%name
else if(present(vfield)) then
call profiler_tic(vfield, "petsc_setup")
name = vfield%name
else if(present(tfield)) then
call profiler_tic(tfield, "petsc_setup")
name = tfield%name
else
FLAbort("petsc_solve_setup should be called with sfield, vfield or tfield")
end if
timing=(debug_level()>=2)
if (timing) then
call cpu_time(time1)
end if
if (present(option_path)) then
solver_option_path=complete_solver_option_path(option_path)
else if (present(sfield)) then
solver_option_path=complete_solver_option_path(sfield%option_path)
else if (present(vfield)) then
solver_option_path=complete_solver_option_path(vfield%option_path)
else if (present(tfield)) then
solver_option_path=complete_solver_option_path(tfield%option_path)
else
FLAbort("Need to provide either sfield, vfield or tfield to petsc_solve_setup.")
end if
startfromzero=have_option(trim(solver_option_path)//'/start_from_zero')
if (present_and_true(startfromzero_in) .and. .not. startfromzero) then
ewrite(2,*) 'Note: startfromzero hard-coded to .true.'
ewrite(2,*) 'Ignoring setting from solver option.'
startfromzero=.true.
end if
ksp=PETSC_NULL_KSP
if (present(matrix)) then
if (associated(matrix%ksp)) then
ksp=matrix%ksp
end if
else if (present(block_matrix)) then
if (associated(block_matrix%ksp)) then
ksp=block_matrix%ksp
end if
end if
if (ksp/=PETSC_NULL_KSP) then
! oh goody, we've been left something useful!
call KSPGetOperators(ksp, A, Pmat, ierr)
have_cache=.true.
if (have_option(trim(solver_option_path)// &
'/preconditioner::mg/vertical_lumping/internal_smoother')) then
! this option is unsafe with caching, as it needs
! DestroyMultigrid to be called on top of PCDestroy to destroy
! all its associated objects
FLExit("Sorry, can't combine internal_smoother with cache_solver_context")
end if
else
! no cache - we just have to do it all over again
have_cache=.false.
end if
ewrite(1, *) 'Assembling matrix.'
! Note the explicitly-described options rcm, 1wd and natural are now not
! listed explicitly in the schema (but can still be used by adding the
! appropriate string in the solver reordering node).
call PetscOptionsGetString(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, "-ordering_type", ordering_type, use_reordering, ierr)
if (.not. use_reordering) then
call get_option(trim(solver_option_path)//'/reordering[0]/name', &
ordering_type, stat=ierr)
use_reordering= (ierr==0)
end if
if (present(matrix)) then
ewrite(2, *) 'Number of rows == ', size(matrix, 1)
! Create the matrix & vectors.
inactive_mask => get_inactive_mask(matrix)
! create list of inactive, ghost_nodes
if (associated(inactive_mask)) then
allocate( ghost_nodes(1:count(inactive_mask)) )
j=0
do i=1, size(matrix,1)
if (inactive_mask(i)) then
j=j+1
ghost_nodes(j)=i
end if
end do
else
allocate( ghost_nodes(1:0) )
end if
! set up numbering used in PETSc objects:
! NOTE: we use size(matrix,2) here as halo rows may be absent
call allocate(petsc_numbering, &
nnodes=size(matrix,2), nfields=1, &
halo=matrix%sparsity%column_halo, &
ghost_nodes=ghost_nodes)
if (use_reordering) then
call reorder(petsc_numbering, matrix%sparsity, ordering_type)
end if
if (.not. have_cache) then
! create PETSc Mat using this numbering:
A=csr2petsc(matrix, petsc_numbering, petsc_numbering)
end if
halo=>matrix%sparsity%column_halo
elseif (present(block_matrix)) then
ewrite(2, *) 'Number of rows == ', size(block_matrix, 1)
ewrite(2, *) 'Number of blocks == ', blocks(block_matrix,1)
assert(.not.block_matrix%diagonal)
! Create the matrix & vectors.
! set up numbering used in PETSc objects:
call allocate(petsc_numbering, &
nnodes=block_size(block_matrix,2), nfields=blocks(block_matrix,1), &
halo=block_matrix%sparsity%column_halo)
if (use_reordering) then
call reorder(petsc_numbering, block_matrix%sparsity, ordering_type)
end if
if (.not. have_cache) then
! create PETSc Mat using this numbering:
A=block_csr2petsc(block_matrix, petsc_numbering, petsc_numbering)
end if
halo=>block_matrix%sparsity%column_halo
else
ewrite(-1,*) "So what am I going to solve???"
FLAbort("Wake up!")
end if
ewrite(1, *) 'Matrix assembly completed.'
if (IsParallel()) then
parallel= (associated(halo))
else
parallel=.false.
end if
if (have_cache) then
! write the cached solver options to log:
call ewrite_ksp_options(ksp)
else
if (present(preconditioner_matrix)) then
ewrite(2,*) 'Using provided preconditioner matrix'
pmat=csr2petsc(preconditioner_matrix, petsc_numbering)
else
pmat=A
end if
ewrite(2, *) 'Using solver options defined at: ', trim(solver_option_path)
call attach_null_space_from_options(A, solver_option_path, pmat=pmat, &
positions=positions, petsc_numbering=petsc_numbering)
call create_ksp_from_options(ksp, A, pmat, solver_option_path, parallel, &
petsc_numbering, &
startfromzero_in=startfromzero_in, &
prolongators=prolongators, surface_node_list=surface_node_list, &
matrix_csr=matrix, &
internal_smoothing_option=internal_smoothing_option)
end if
if (.not. have_cache .and. have_option(trim(solver_option_path)// &
&'/cache_solver_context')) then
! save the ksp solver context for future generations
! (hack with pointer to convince intel compiler that it's
! really just the pointed-to value I'm changing)
if (present(matrix)) then
ksp_pointer => matrix%ksp
else if (present(block_matrix)) then
ksp_pointer => block_matrix%ksp
end if
if (associated(ksp_pointer)) then
ksp_pointer = ksp
! make sure we don't destroy it, the %ksp becomes a separate reference
call PetscObjectReferenceWrapper(ksp, ierr)
else
! matrices coming from block() can't cache
FLAbort("User wants to cache solver context, but no proper matrix is provided.")
end if
else if (have_cache) then
! ksp is a copy of matrix%ksp, make it a separate reference,
! so we can KSPDestroy it without destroying matrix%ksp
call PetscObjectReferenceWrapper(ksp, ierr)
! same for the matrix, kspgetoperators returns the matrix reference
! owned by the ksp - make it a separate reference
call PetscObjectReferenceWrapper(A, ierr)
end if
b=PetscNumberingCreateVec(petsc_numbering)
call FixedVecDuplicate(b, y, ierr)
if (timing) then
call cpu_time(time2)
ewrite(2,*) trim(name)// " CPU time spent in PETSc setup: ", time2-time1
end if
if(present(sfield)) then
call profiler_toc(sfield, "petsc_setup")
else if(present(vfield)) then
call profiler_toc(vfield, "petsc_setup")
else if(present(tfield)) then
call profiler_toc(tfield, "petsc_setup")
end if
end subroutine petsc_solve_setup
subroutine petsc_solve_setup_petsc_csr(y, b, &
solver_option_path, startfromzero, &
matrix, sfield, vfield, tfield, &
option_path, startfromzero_in, &
prolongators,surface_node_list, &
positions, rotation_matrix)
!!< sets up things needed to call petsc_solve_core
!! Stuff that comes out:
!!
!! PETSc solution vector
Vec, intent(out):: y
!! PETSc rhs vector
Vec, intent(out):: b
!! returns the option path to solver/ block for new options, otherwise ""
character(len=*), intent(out):: solver_option_path
!! whether to start with zero initial guess
logical, intent(out):: startfromzero
!! Stuff that goes in:
!!
!! provide either a matrix or block_matrix to be solved
type(petsc_csr_matrix), intent(inout):: matrix
!! provide either a scalar field or vector field to be solved for
type(scalar_field), optional, intent(in):: sfield
type(vector_field), optional, intent(in):: vfield
type(tensor_field), optional, intent(in):: tfield
!! overrides sfield%option_path or vfield%option_path
character(len=*), intent(in), optional:: option_path
!! whether to start with zero initial guess (as passed in)
logical, optional, intent(in):: startfromzero_in
!! additional info for "mg" preconditioner:
!! prolongators to be used at the first levels of 'mg'
type(petsc_csr_matrix), dimension(:), optional, intent(in) :: prolongators
!! Stuff needed for internal smoother
integer, dimension(:), optional, intent(in) :: surface_node_list
!! positions field is only used with remove_null_space/ or multigrid_near_null_space/ with rotational components
type(vector_field), intent(in), optional :: positions
!! with rotated bcs: matrix to transform from x,y,z aligned vectors to boundary aligned
Mat, intent(in), optional:: rotation_matrix
real time1, time2
integer ierr
logical parallel, timing
character(len=FIELD_NAME_LEN) :: name
! Initialise profiler
if(present(sfield)) then
call profiler_tic(sfield, "petsc_setup")
name = sfield%name
else if(present(vfield)) then
call profiler_tic(vfield, "petsc_setup")
name = vfield%name
else if(present(tfield)) then
call profiler_tic(tfield, "petsc_setup")
name = tfield%name
else
FLAbort("petsc_solve_setup should be called with sfield, vfield or tfield")
end if
timing=(debug_level()>=2)
if (timing) then
call cpu_time(time1)
end if
if (present(option_path)) then
solver_option_path=complete_solver_option_path(option_path)
else if (present(sfield)) then
solver_option_path=complete_solver_option_path(sfield%option_path)
else if (present(vfield)) then
solver_option_path=complete_solver_option_path(vfield%option_path)
else if (present(tfield)) then
solver_option_path=complete_solver_option_path(tfield%option_path)
else
FLAbort("Need to provide either sfield, vfield or tfield to petsc_solve_setup.")
end if