File: kernel.html

package info (click to toggle)
petsc 3.23.1%2Bdfsg1-1exp1
  • links: PTS, VCS
  • area: main
  • in suites: experimental
  • size: 515,576 kB
  • sloc: ansic: 751,607; cpp: 51,542; python: 38,598; f90: 17,352; javascript: 3,493; makefile: 3,157; sh: 1,502; xml: 619; objc: 445; java: 13; csh: 1
file content (1157 lines) | stat: -rw-r--r-- 118,070 bytes parent folder | download
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
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157

<!DOCTYPE html>


<html lang="en" data-content_root="../" >

  <head>
    <meta charset="utf-8" />
    <meta name="viewport" content="width=device-width, initial-scale=1.0" /><meta name="viewport" content="width=device-width, initial-scale=1" />

    <title>The PETSc Kernel &#8212; PETSc 3.23.1 documentation</title>
  
  
  
  <script data-cfasync="false">
    document.documentElement.dataset.mode = localStorage.getItem("mode") || "";
    document.documentElement.dataset.theme = localStorage.getItem("theme") || "light";
  </script>
  
  <!-- Loaded before other Sphinx assets -->
  <link href="../_static/styles/theme.css?digest=bd9e20870c6007c4c509" rel="stylesheet" />
<link href="../_static/styles/bootstrap.css?digest=bd9e20870c6007c4c509" rel="stylesheet" />
<link href="../_static/styles/pydata-sphinx-theme.css?digest=bd9e20870c6007c4c509" rel="stylesheet" />

  
  <link href="../_static/vendor/fontawesome/6.5.1/css/all.min.css?digest=bd9e20870c6007c4c509" rel="stylesheet" />
  <link rel="preload" as="font" type="font/woff2" crossorigin href="../_static/vendor/fontawesome/6.5.1/webfonts/fa-solid-900.woff2" />
<link rel="preload" as="font" type="font/woff2" crossorigin href="../_static/vendor/fontawesome/6.5.1/webfonts/fa-brands-400.woff2" />
<link rel="preload" as="font" type="font/woff2" crossorigin href="../_static/vendor/fontawesome/6.5.1/webfonts/fa-regular-400.woff2" />

    <link rel="stylesheet" type="text/css" href="../_static/pygments.css?v=8f2a1f02" />
    <link rel="stylesheet" type="text/css" href="../_static/copybutton.css?v=76b2166b" />
    <link rel="stylesheet" type="text/css" href="../_static/sphinx-design.min.css?v=87e54e7c" />
    <link rel="stylesheet" type="text/css" href="https://cdn.jsdelivr.net/npm/katex@0.16.10/dist/katex.min.css" />
    <link rel="stylesheet" type="text/css" href="../_static/katex-math.css?v=91adb8b6" />
    <link rel="stylesheet" type="text/css" href="../_static/css/custom.css?v=dbe1606d" />
  
  <!-- Pre-loaded scripts that we'll load fully later -->
  <link rel="preload" as="script" href="../_static/scripts/bootstrap.js?digest=bd9e20870c6007c4c509" />
<link rel="preload" as="script" href="../_static/scripts/pydata-sphinx-theme.js?digest=bd9e20870c6007c4c509" />
  <script src="../_static/vendor/fontawesome/6.5.1/js/all.min.js?digest=bd9e20870c6007c4c509"></script>

    <script src="../_static/documentation_options.js?v=34da53a5"></script>
    <script src="../_static/doctools.js?v=9a2dae69"></script>
    <script src="../_static/sphinx_highlight.js?v=dc90522c"></script>
    <script src="../_static/clipboard.min.js?v=a7894cd8"></script>
    <script src="../_static/copybutton.js?v=a56c686a"></script>
    <script src="../_static/design-tabs.js?v=f930bc37"></script>
    <script src="../_static/katex.min.js?v=be8ff15f"></script>
    <script src="../_static/auto-render.min.js?v=ad136472"></script>
    <script src="../_static/katex_autorenderer.js?v=bebc588a"></script>
    <script>DOCUMENTATION_OPTIONS.pagename = 'developers/kernel';</script>
    <link rel="icon" href="../_static/petsc_favicon.png"/>
    <link rel="index" title="Index" href="../genindex.html" />
    <link rel="search" title="Search" href="../search.html" />
    <link rel="next" title="Basic Object Design and Implementation" href="objects.html" />
    <link rel="prev" title="The Design of PETSc" href="design.html" />
  <meta name="viewport" content="width=device-width, initial-scale=1"/>
  <meta name="docsearch:language" content="en"/>
    <meta name="docbuild:last-update" content="2025-04-30T13:10:40-0500 (v3.23.1)"/>
  </head>
  
  
  <body data-bs-spy="scroll" data-bs-target=".bd-toc-nav" data-offset="180" data-bs-root-margin="0px 0px -60%" data-default-mode="">

  
  
  <a id="pst-skip-link" class="skip-link" href="#main-content">Skip to main content</a>
  
  <div id="pst-scroll-pixel-helper"></div>

  
  <button type="button" class="btn rounded-pill" id="pst-back-to-top">
    <i class="fa-solid fa-arrow-up"></i>
    Back to top
  </button>

  
  <input type="checkbox"
          class="sidebar-toggle"
          name="__primary"
          id="__primary"/>
  <label class="overlay overlay-primary" for="__primary"></label>
  
  <input type="checkbox"
          class="sidebar-toggle"
          name="__secondary"
          id="__secondary"/>
  <label class="overlay overlay-secondary" for="__secondary"></label>
  
  <div class="search-button__wrapper">
    <div class="search-button__overlay"></div>
    <div class="search-button__search-container">
<form class="bd-search d-flex align-items-center"
      action="../search.html"
      method="get">
  <i class="fa-solid fa-magnifying-glass"></i>
  <input type="search"
         class="form-control"
         name="q"
         id="search-input"
         placeholder="Search the docs ..."
         aria-label="Search the docs ..."
         autocomplete="off"
         autocorrect="off"
         autocapitalize="off"
         spellcheck="false"/>
  <span class="search-button__kbd-shortcut"><kbd class="kbd-shortcut__modifier">Ctrl</kbd>+<kbd>K</kbd></span>
</form></div>
  </div>

  <header>
  
    <div class="bd-header navbar navbar-expand-lg bd-navbar">
<div class="bd-header__inner bd-page-width">
  <label class="sidebar-toggle primary-toggle" for="__primary">
    <span class="fa-solid fa-bars"></span>
  </label>
  
  
  <div class="col-lg-3 navbar-header-items__start">
    
      <div class="navbar-item">

  

<a class="navbar-brand logo" href="../index.html">
  
  
  
  
  
    
    
      
    
    
    <img src="../_static/PETSc-TAO_RGB.svg" class="logo__image only-light" alt="PETSc 3.23.1 documentation - Home"/>
    <script>document.write(`<img src="../_static/PETSc-TAO_RGB_white.svg" class="logo__image only-dark" alt="PETSc 3.23.1 documentation - Home"/>`);</script>
  
  
</a></div>
    
  </div>
  
  <div class="col-lg-9 navbar-header-items">
    
    <div class="me-auto navbar-header-items__center">
      
        <div class="navbar-item">
<nav class="navbar-nav">
  <ul class="bd-navbar-elements navbar-nav">
    
                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../overview/index.html">
                        Overview
                      </a>
                    </li>
                

                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../install/index.html">
                        Install
                      </a>
                    </li>
                

                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../tutorials/index.html">
                        Tutorials
                      </a>
                    </li>
                

                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../manual/index.html">
                        User-Guide
                      </a>
                    </li>
                

                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../manualpages/index.html">
                        C/Fortran API
                      </a>
                    </li>
                

                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../petsc4py/index.html">
                        petsc4py API
                      </a>
                    </li>
                

                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../faq/index.html">
                        FAQ
                      </a>
                    </li>
                

                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../community/index.html">
                        Community
                      </a>
                    </li>
                

                    <li class="nav-item current active">
                      <a class="nav-link nav-internal" href="index.html">
                        Developers
                      </a>
                    </li>
                

                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../miscellaneous/index.html">
                        Misc.
                      </a>
                    </li>
                
  </ul>
</nav></div>
      
    </div>
    
    
    <div class="navbar-header-items__end">
      
        <div class="navbar-item navbar-persistent--container">
          

 <script>
 document.write(`
   <button class="btn navbar-btn search-button-field search-button__button" title="Search" aria-label="Search" data-bs-placement="bottom" data-bs-toggle="tooltip">
    <i class="fa-solid fa-magnifying-glass"></i>
    <span class="search-button__default-text">Search</span>
    <span class="search-button__kbd-shortcut"><kbd class="kbd-shortcut__modifier">Ctrl</kbd>+<kbd class="kbd-shortcut__modifier">K</kbd></span>
   </button>
 `);
 </script>
        </div>
      
      
        <div class="navbar-item">

<script>
document.write(`
  <button class="btn btn-sm navbar-btn theme-switch-button" title="light/dark" aria-label="light/dark" data-bs-placement="bottom" data-bs-toggle="tooltip">
    <span class="theme-switch nav-link" data-mode="light"><i class="fa-solid fa-sun fa-lg"></i></span>
    <span class="theme-switch nav-link" data-mode="dark"><i class="fa-solid fa-moon fa-lg"></i></span>
    <span class="theme-switch nav-link" data-mode="auto"><i class="fa-solid fa-circle-half-stroke fa-lg"></i></span>
  </button>
`);
</script></div>
      
        <div class="navbar-item"><ul class="navbar-icon-links navbar-nav"
    aria-label="Icon Links">
        <li class="nav-item">
          
          
          
          
          
          
          
          
          <a href="https://gitlab.com/petsc/petsc" title="GitLab" class="nav-link" rel="noopener" target="_blank" data-bs-toggle="tooltip" data-bs-placement="bottom"><span><i class="fab fa-gitlab fa-lg" aria-hidden="true"></i></span>
            <span class="sr-only">GitLab</span></a>
        </li>
</ul></div>
      
    </div>
    
  </div>
  
  
    <div class="navbar-persistent--mobile">

 <script>
 document.write(`
   <button class="btn navbar-btn search-button-field search-button__button" title="Search" aria-label="Search" data-bs-placement="bottom" data-bs-toggle="tooltip">
    <i class="fa-solid fa-magnifying-glass"></i>
    <span class="search-button__default-text">Search</span>
    <span class="search-button__kbd-shortcut"><kbd class="kbd-shortcut__modifier">Ctrl</kbd>+<kbd class="kbd-shortcut__modifier">K</kbd></span>
   </button>
 `);
 </script>
    </div>
  

  
    <label class="sidebar-toggle secondary-toggle" for="__secondary" tabindex="0">
      <span class="fa-solid fa-outdent"></span>
    </label>
  
</div>

    </div>
  
  </header>

  <div class="bd-container">
    <div class="bd-container__inner bd-page-width">
      
      
      
      <div class="bd-sidebar-primary bd-sidebar">
        

  
  <div class="sidebar-header-items sidebar-primary__section">
    
    
      <div class="sidebar-header-items__center">
        
          <div class="navbar-item">
<nav class="navbar-nav">
  <ul class="bd-navbar-elements navbar-nav">
    
                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../overview/index.html">
                        Overview
                      </a>
                    </li>
                

                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../install/index.html">
                        Install
                      </a>
                    </li>
                

                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../tutorials/index.html">
                        Tutorials
                      </a>
                    </li>
                

                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../manual/index.html">
                        User-Guide
                      </a>
                    </li>
                

                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../manualpages/index.html">
                        C/Fortran API
                      </a>
                    </li>
                

                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../petsc4py/index.html">
                        petsc4py API
                      </a>
                    </li>
                

                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../faq/index.html">
                        FAQ
                      </a>
                    </li>
                

                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../community/index.html">
                        Community
                      </a>
                    </li>
                

                    <li class="nav-item current active">
                      <a class="nav-link nav-internal" href="index.html">
                        Developers
                      </a>
                    </li>
                

                    <li class="nav-item">
                      <a class="nav-link nav-internal" href="../miscellaneous/index.html">
                        Misc.
                      </a>
                    </li>
                
  </ul>
</nav></div>
        
      </div>
    
    
    
      <div class="sidebar-header-items__end">
        
          <div class="navbar-item">

<script>
document.write(`
  <button class="btn btn-sm navbar-btn theme-switch-button" title="light/dark" aria-label="light/dark" data-bs-placement="bottom" data-bs-toggle="tooltip">
    <span class="theme-switch nav-link" data-mode="light"><i class="fa-solid fa-sun fa-lg"></i></span>
    <span class="theme-switch nav-link" data-mode="dark"><i class="fa-solid fa-moon fa-lg"></i></span>
    <span class="theme-switch nav-link" data-mode="auto"><i class="fa-solid fa-circle-half-stroke fa-lg"></i></span>
  </button>
`);
</script></div>
        
          <div class="navbar-item"><ul class="navbar-icon-links navbar-nav"
    aria-label="Icon Links">
        <li class="nav-item">
          
          
          
          
          
          
          
          
          <a href="https://gitlab.com/petsc/petsc" title="GitLab" class="nav-link" rel="noopener" target="_blank" data-bs-toggle="tooltip" data-bs-placement="bottom"><span><i class="fab fa-gitlab fa-lg" aria-hidden="true"></i></span>
            <span class="sr-only">GitLab</span></a>
        </li>
</ul></div>
        
      </div>
    
  </div>
  
    <div class="sidebar-primary-items__start sidebar-primary__section">
        <div class="sidebar-primary-item">
<nav class="bd-docs-nav bd-links"
     aria-label="Section Navigation">
  <p class="bd-links__title" role="heading" aria-level="1">Section Navigation</p>
  <div class="bd-toc-item navbar-nav"><ul class="current nav bd-sidenav">
<li class="toctree-l1"><a class="reference internal" href="communication.html">PETSc Developers Communication Channels</a></li>
<li class="toctree-l1 has-children"><a class="reference internal" href="contributing/index.html">Contributing to PETSc</a><input class="toctree-checkbox" id="toctree-checkbox-1" name="toctree-checkbox-1" type="checkbox"/><label class="toctree-toggle" for="toctree-checkbox-1"><i class="fa-solid fa-chevron-down"></i></label><ul>
<li class="toctree-l2"><a class="reference internal" href="contributing/developingmr.html">Developing a Merge Request</a></li>
<li class="toctree-l2"><a class="reference internal" href="contributing/submittingmr.html">Submitting a Merge Request</a></li>
<li class="toctree-l2"><a class="reference internal" href="contributing/pipelines.html">GitLab CI Pipelines</a></li>
</ul>
</li>
<li class="toctree-l1"><a class="reference internal" href="mrmanagement.html">Merge request management</a></li>
<li class="toctree-l1"><a class="reference internal" href="development.html">PETSc Development Environment</a></li>
<li class="toctree-l1"><a class="reference internal" href="style.html">PETSc Style and Usage Guide</a></li>
<li class="toctree-l1"><a class="reference internal" href="buildsystem.html">BuildSystem</a></li>
<li class="toctree-l1"><a class="reference internal" href="testing.html">PETSc Testing System</a></li>
<li class="toctree-l1"><a class="reference internal" href="documentation.html">Developing PETSc Documentation</a></li>
<li class="toctree-l1 current active has-children"><a class="reference internal" href="design.html">The Design of PETSc</a><input checked="" class="toctree-checkbox" id="toctree-checkbox-2" name="toctree-checkbox-2" type="checkbox"/><label class="toctree-toggle" for="toctree-checkbox-2"><i class="fa-solid fa-chevron-down"></i></label><ul class="current">
<li class="toctree-l2 current active"><a class="current reference internal" href="#">The PETSc Kernel</a></li>
<li class="toctree-l2"><a class="reference internal" href="objects.html">Basic Object Design and Implementation</a></li>
<li class="toctree-l2"><a class="reference internal" href="callbacks.html">How the Solvers Handle User Provided Callbacks</a></li>
<li class="toctree-l2"><a class="reference internal" href="matrices.html">The Various Matrix Classes</a></li>
<li class="toctree-l2"><a class="reference internal" href="articles.html">Articles about PETSc Design</a></li>
</ul>
</li>
</ul>
</div>
</nav></div>
    </div>
  
  
  <div class="sidebar-primary-items__end sidebar-primary__section">
  </div>
  
  <div id="rtd-footer-container"></div>


      </div>
      
      <main id="main-content" class="bd-main">
        
        
          <div class="bd-content">
            <div class="bd-article-container">
              
              <div class="bd-header-article">
<div class="header-article-items header-article__inner">
  
    <div class="header-article-items__start">
      
        <div class="header-article-item">



<nav aria-label="Breadcrumb">
  <ul class="bd-breadcrumbs">
    
    <li class="breadcrumb-item breadcrumb-home">
      <a href="../index.html" class="nav-link" aria-label="Home">
        <i class="fa-solid fa-home"></i>
      </a>
    </li>
    
    <li class="breadcrumb-item"><a href="index.html" class="nav-link">Developers</a></li>
    
    
    <li class="breadcrumb-item"><a href="design.html" class="nav-link">The Design of PETSc</a></li>
    
    <li class="breadcrumb-item active" aria-current="page">The PETSc Kernel</li>
  </ul>
</nav>
</div>
      
    </div>
  
  
</div>
</div>
              
              
              
                
<div id="searchbox"></div>
                <article class="bd-article">
                  
  <section class="tex2jax_ignore mathjax_ignore" id="the-petsc-kernel">
<h1>The PETSc Kernel<a class="headerlink" href="#the-petsc-kernel" title="Link to this heading">#</a></h1>
<p>PETSc provides a variety of basic services for writing scalable,
component-based libraries; these are referred to as the PETSc kernel <span id="id1">[<a class="reference internal" href="#id821" title="Satish Balay, William D. Gropp, Lois Curfman McInnes, and Barry F. Smith. A microkernel design for component-based parallel numerical software systems. In Proceedings of the SIAM Workshop on Object Oriented Methods for Inter-operable Scientific and Engineering Computing, 58-67. SIAM, 1998. URL: ftp://info.mcs.anl.gov/pub/tech_reports/reports/P727.ps.Z.">BGMS98</a>]</span>.
The source code for the kernel is in <code class="docutils notranslate"><span class="pre">src/sys</span></code>. It contains systematic
support for</p>
<ul class="simple">
<li><p>managing PETSc types,</p></li>
<li><p>error handling,</p></li>
<li><p>memory management,</p></li>
<li><p>profiling,</p></li>
<li><p>object management,</p></li>
<li><p>Fortran interfaces (see <span id="id2">[<a class="reference internal" href="#id958" title="Satish Balay, Jed Brown, Matthew G. Knepley, Lois McInnes, and Barry Smith. Providing mixed language and legacy support within a library. In J. Carver, editor, Software Engineering for Science. Taylor &amp; Francis, 2015.">BBK+15</a>]</span>)</p></li>
<li><p>mechanism for generating appropriate citations for algorithms and software used in PETSc (see <span id="id3">[<a class="reference internal" href="#id960" title="Matthew G. Knepley, Jed Brown, Lois Curfman McInnes, and Barry F. Smith. Accurately citing software and algorithms used in publications. In First Workshop on Sustainable Software for Science: Practice and Experiences (WSSSPE), held at SC13. 2013. doi:10.6084/m9.figshare.785731.">KBMS13</a>]</span>)</p></li>
<li><p>file I/O,</p></li>
<li><p>an options database, and</p></li>
<li><p>objects and code for viewing, drawing, and displaying data and solver objects.</p></li>
</ul>
<p>Each of these is discussed in a section below.</p>
<section id="petsc-types">
<h2>PETSc Types<a class="headerlink" href="#petsc-types" title="Link to this heading">#</a></h2>
<p>For maximum flexibility, the basic data types <code class="docutils notranslate"><span class="pre">int</span></code>, <code class="docutils notranslate"><span class="pre">double</span></code>, and
so on are not used in PETSc source code. Rather, it has</p>
<ul class="simple">
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a></span></code>,</p></li>
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscInt.html">PetscInt</a></span></code>,</p></li>
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscMPIInt.html">PetscMPIInt</a></span></code>,</p></li>
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscBLASInt.html">PetscBLASInt</a></span></code>,</p></li>
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscBool.html">PetscBool</a></span></code>, and</p></li>
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscBT.html">PetscBT</a></span></code> - bit storage of logical true and false.</p></li>
</ul>
<p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscInt.html">PetscInt</a></span></code> can be set using <code class="docutils notranslate"><span class="pre">configure</span></code> to be either <code class="docutils notranslate"><span class="pre">int</span></code> (32
bit, the default) or <code class="docutils notranslate"><span class="pre">long</span> <span class="pre">long</span></code> (64-bit, with
<code class="docutils notranslate"><span class="pre">configure</span> <span class="pre">–with-64-bit-indices</span></code>) to allow indexing into very large
arrays. <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscMPIInt.html">PetscMPIInt</a></span></code> is used for integers passed to MPI as counts and
sizes. These are always <code class="docutils notranslate"><span class="pre">int</span></code> since that is what the MPI standard
uses. Similarly, <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscBLASInt.html">PetscBLASInt</a></span></code> is for counts, and so on passed to
BLAS and LAPACK routines. These are almost always <code class="docutils notranslate"><span class="pre">int</span></code> unless one is
using a special “64-bit integer” BLAS/LAPACK (this is available, for
example, with Intel’s MKL and OpenBLAS).</p>
<p>In addition, there are special types:</p>
<ul class="simple">
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscClassId.html">PetscClassId</a></span></code></p></li>
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a></span></code></p></li>
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Log/PetscLogEvent.html">PetscLogEvent</a></span></code></p></li>
</ul>
<p>These are currently always <code class="docutils notranslate"><span class="pre">int</span></code>, but their use clarifies the code.</p>
</section>
<section id="implementation-of-error-handling">
<h2>Implementation of Error Handling<a class="headerlink" href="#implementation-of-error-handling" title="Link to this heading">#</a></h2>
<p>PETSc uses a “call error handler; then (depending on result) return
error code” model when problems are detected in the running code. The
public include file for error handling is
<a href="../include/petscerror.h.html">include/petscerror.h</a>
and the source code for the PETSc error handling is in
<code class="docutils notranslate"><span class="pre">src/sys/error/</span></code>.</p>
<section id="simplified-interface">
<h3>Simplified Interface<a class="headerlink" href="#simplified-interface" title="Link to this heading">#</a></h3>
<p>The simplified macro-based interface consists of the following two
calls:</p>
<ul class="simple">
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/SETERRQ.html">SETERRQ</a>(comm,error</span> <span class="pre">code,Error</span> <span class="pre">message);</span></code></p></li>
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a>(ierr);</span></code></p></li>
</ul>
<p>The macro <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/SETERRQ.html">SETERRQ</a>()</span></code> is given by</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="cp">#define SETERRQ1(...) PETSC_DEPRECATED_MACRO(3, 17, 0, &quot;<a href="../manualpages/Sys/SETERRQ.html">SETERRQ</a>&quot;, ) <a href="../manualpages/Sys/SETERRQ.html">SETERRQ</a>(__VA_ARGS__)</span>
</pre></div>
</div>
<p>It calls the error handler with the current function name and location:
line number, and file, plus an error code and an error message. Normally
<code class="docutils notranslate"><span class="pre">comm</span></code> is <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PETSC_COMM_SELF.html">PETSC_COMM_SELF</a></span></code>; it can be another communicator only if
one is absolutely sure the same error will be generated on all processes
in the communicator. This feature is to prevent the same error message
from being printed by many processes.</p>
<p>The macro <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a>()</span></code> is defined by</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="w">  </span><span class="cp">#define <a href="../manualpages/Sys/PetscCall.html">PetscCall</a>(...) \</span>
</pre></div>
</div>
<p>The message passed to <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/SETERRQ.html">SETERRQ</a>()</span></code> is treated as a <code class="docutils notranslate"><span class="pre">printf()</span></code>-style
format string, with all additional parameters passed after the string as
its arguments. For example:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../manualpages/Sys/SETERRQ.html">SETERRQ</a></span><span class="p">(</span><span class="n">comm</span><span class="p">,</span><span class="n">PETSC_ERR</span><span class="p">,</span><span class="s">&quot;Iteration overflow: its %&quot;</span><span class="w"> </span><span class="n">PetscInt_FMT</span><span class="w"> </span><span class="s">&quot; norm %g&quot;</span><span class="p">,</span><span class="n">its</span><span class="p">,(</span><span class="kt">double</span><span class="p">)</span><span class="n">norm</span><span class="p">);</span>
</pre></div>
</div>
</section>
<section id="error-handlers">
<h3>Error Handlers<a class="headerlink" href="#error-handlers" title="Link to this heading">#</a></h3>
<p>The error-handling function <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscError.html">PetscError</a>()</span></code> calls the “current” error
handler with the code</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a></span><span class="w"> </span><span class="nf"><a href="../manualpages/Sys/PetscError.html">PetscError</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/MPI_Comm.html">MPI_Comm</a></span><span class="w"> </span><span class="n">comm</span><span class="p">,</span><span class="w"> </span><span class="kt">int</span><span class="w"> </span><span class="n">line</span><span class="p">,</span><span class="w"> </span><span class="k">const</span><span class="w"> </span><span class="kt">char</span><span class="w"> </span><span class="o">*</span><span class="n">func</span><span class="p">,</span><span class="w"> </span><span class="k">const</span><span class="w"> </span><span class="kt">char</span><span class="w"> </span><span class="o">*</span><span class="n">file</span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a></span><span class="w"> </span><span class="n">n</span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscErrorType.html">PetscErrorType</a></span><span class="w"> </span><span class="n">p</span><span class="p">,</span><span class="w"> </span><span class="k">const</span><span class="w"> </span><span class="kt">char</span><span class="w"> </span><span class="o">*</span><span class="n">mess</span><span class="p">,</span><span class="w"> </span><span class="p">...)</span>
<span class="p">{</span>
<span class="w">  </span><span class="kt">va_list</span><span class="w">        </span><span class="n">Argp</span><span class="p">;</span>
<span class="w">  </span><span class="kt">size_t</span><span class="w">         </span><span class="n">fullLength</span><span class="p">;</span>
<span class="w">  </span><span class="kt">char</span><span class="w">           </span><span class="n">buf</span><span class="p">[</span><span class="mi">2048</span><span class="p">],</span><span class="w"> </span><span class="o">*</span><span class="n">lbuf</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="nb">NULL</span><span class="p">;</span>
<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscBool.html">PetscBool</a></span><span class="w">      </span><span class="n">ismain</span><span class="p">;</span>
<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a></span><span class="w"> </span><span class="n">ierr</span><span class="p">;</span>

<span class="w">  </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="o">!</span><span class="n">PetscErrorHandlingInitialized</span><span class="p">)</span><span class="w"> </span><span class="k">return</span><span class="w"> </span><span class="n">n</span><span class="p">;</span>
<span class="w">  </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">comm</span><span class="w"> </span><span class="o">==</span><span class="w"> </span><span class="n">MPI_COMM_NULL</span><span class="p">)</span><span class="w"> </span><span class="n">comm</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PETSC_COMM_SELF.html">PETSC_COMM_SELF</a></span><span class="p">;</span>

<span class="w">  </span><span class="cm">/* Compose the message evaluating the print format */</span>
<span class="w">  </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">mess</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">    </span><span class="n">va_start</span><span class="p">(</span><span class="n">Argp</span><span class="p">,</span><span class="w"> </span><span class="n">mess</span><span class="p">);</span>
<span class="w">    </span><span class="p">(</span><span class="kt">void</span><span class="p">)</span><span class="n"><a href="../manualpages/Sys/PetscVSNPrintf.html">PetscVSNPrintf</a></span><span class="p">(</span><span class="n">buf</span><span class="p">,</span><span class="w"> </span><span class="mi">2048</span><span class="p">,</span><span class="w"> </span><span class="n">mess</span><span class="p">,</span><span class="w"> </span><span class="o">&amp;</span><span class="n">fullLength</span><span class="p">,</span><span class="w"> </span><span class="n">Argp</span><span class="p">);</span>
<span class="w">    </span><span class="n">va_end</span><span class="p">(</span><span class="n">Argp</span><span class="p">);</span>
<span class="w">    </span><span class="n">lbuf</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="n">buf</span><span class="p">;</span>
<span class="w">    </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">p</span><span class="w"> </span><span class="o">==</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscErrorType.html">PETSC_ERROR_INITIAL</a></span><span class="p">)</span><span class="w"> </span><span class="p">(</span><span class="kt">void</span><span class="p">)</span><span class="n"><a href="../manualpages/Sys/PetscStrncpy.html">PetscStrncpy</a></span><span class="p">(</span><span class="n">PetscErrorBaseMessage</span><span class="p">,</span><span class="w"> </span><span class="n">lbuf</span><span class="p">,</span><span class="w"> </span><span class="k">sizeof</span><span class="p">(</span><span class="n">PetscErrorBaseMessage</span><span class="p">));</span>
<span class="w">  </span><span class="p">}</span>

<span class="w">  </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">p</span><span class="w"> </span><span class="o">==</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscErrorType.html">PETSC_ERROR_INITIAL</a></span><span class="w"> </span><span class="o">&amp;&amp;</span><span class="w"> </span><span class="n">n</span><span class="w"> </span><span class="o">!=</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscErrorCode.html">PETSC_ERR_MEMC</a></span><span class="p">)</span><span class="w"> </span><span class="p">(</span><span class="kt">void</span><span class="p">)</span><span class="n"><a href="../manualpages/Sys/PetscMallocValidate.html">PetscMallocValidate</a></span><span class="p">(</span><span class="n">__LINE__</span><span class="p">,</span><span class="w"> </span><span class="n">PETSC_FUNCTION_NAME</span><span class="p">,</span><span class="w"> </span><span class="n">__FILE__</span><span class="p">);</span>

<span class="w">  </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="o">!</span><span class="n">eh</span><span class="p">)</span><span class="w"> </span><span class="n">ierr</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscTraceBackErrorHandler.html">PetscTraceBackErrorHandler</a></span><span class="p">(</span><span class="n">comm</span><span class="p">,</span><span class="w"> </span><span class="n">line</span><span class="p">,</span><span class="w"> </span><span class="n">func</span><span class="p">,</span><span class="w"> </span><span class="n">file</span><span class="p">,</span><span class="w"> </span><span class="n">n</span><span class="p">,</span><span class="w"> </span><span class="n">p</span><span class="p">,</span><span class="w"> </span><span class="n">lbuf</span><span class="p">,</span><span class="w"> </span><span class="nb">NULL</span><span class="p">);</span>
<span class="w">  </span><span class="k">else</span><span class="w"> </span><span class="n">ierr</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="p">(</span><span class="o">*</span><span class="n">eh</span><span class="o">-&gt;</span><span class="n">handler</span><span class="p">)(</span><span class="n">comm</span><span class="p">,</span><span class="w"> </span><span class="n">line</span><span class="p">,</span><span class="w"> </span><span class="n">func</span><span class="p">,</span><span class="w"> </span><span class="n">file</span><span class="p">,</span><span class="w"> </span><span class="n">n</span><span class="p">,</span><span class="w"> </span><span class="n">p</span><span class="p">,</span><span class="w"> </span><span class="n">lbuf</span><span class="p">,</span><span class="w"> </span><span class="n">eh</span><span class="o">-&gt;</span><span class="n">ctx</span><span class="p">);</span>
<span class="w">  </span><span class="n">PetscStackClearTop</span><span class="p">;</span>

<span class="w">  </span><span class="cm">/*</span>
<span class="cm">      If this is called from the main() routine we abort the program.</span>
<span class="cm">      We cannot just return because them some MPI processes may continue to attempt to run</span>
<span class="cm">      while this process simply exits.</span>
<span class="cm">  */</span>
<span class="w">  </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">func</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">    </span><span class="p">(</span><span class="kt">void</span><span class="p">)</span><span class="n"><a href="../manualpages/Sys/PetscStrncmp.html">PetscStrncmp</a></span><span class="p">(</span><span class="n">func</span><span class="p">,</span><span class="w"> </span><span class="s">&quot;main&quot;</span><span class="p">,</span><span class="w"> </span><span class="mi">4</span><span class="p">,</span><span class="w"> </span><span class="o">&amp;</span><span class="n">ismain</span><span class="p">);</span>
<span class="w">    </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">ismain</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">      </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">petscwaitonerrorflg</span><span class="p">)</span><span class="w"> </span><span class="p">(</span><span class="kt">void</span><span class="p">)</span><span class="n"><a href="../manualpages/Sys/PetscSleep.html">PetscSleep</a></span><span class="p">(</span><span class="mi">1000</span><span class="p">);</span>
<span class="w">      </span><span class="n"><a href="../manualpages/Sys/PETSCABORT.html">PETSCABORT</a></span><span class="p">(</span><span class="n">comm</span><span class="p">,</span><span class="w"> </span><span class="n">ierr</span><span class="p">);</span>
<span class="w">    </span><span class="p">}</span>
<span class="w">  </span><span class="p">}</span>
<span class="cp">#if defined(PETSC_CLANGUAGE_CXX)</span>
<span class="w">  </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">p</span><span class="w"> </span><span class="o">==</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscErrorType.html">PETSC_ERROR_IN_CXX</a></span><span class="p">)</span><span class="w"> </span><span class="n">PetscCxxErrorThrow</span><span class="p">();</span>
<span class="cp">#endif</span>
<span class="w">  </span><span class="k">return</span><span class="w"> </span><span class="n">ierr</span><span class="p">;</span>
<span class="p">}</span>

<span class="cm">/*@</span>
<span class="cm">  <a href="../manualpages/Sys/PetscIntViewNumColumns.html">PetscIntViewNumColumns</a> - Prints an array of integers; useful for debugging.</span>

<span class="cm">  Collective</span>

<span class="cm">  Input Parameters:</span>
<span class="cm">+ N      - number of integers in array</span>
<span class="cm">. Ncol   - number of integers to print per row</span>
<span class="cm">. idx    - array of integers</span>
<span class="cm">- viewer - location to print array, `<a href="../manualpages/Viewer/PETSC_VIEWER_STDOUT_WORLD.html">PETSC_VIEWER_STDOUT_WORLD</a>`, `<a href="../manualpages/Viewer/PETSC_VIEWER_STDOUT_SELF.html">PETSC_VIEWER_STDOUT_SELF</a>` or 0</span>

<span class="cm">  Level: intermediate</span>

<span class="cm">  Note:</span>
<span class="cm">  This may be called from within the debugger, passing 0 as the viewer</span>

<span class="cm">  This API may be removed in the future.</span>

<span class="cm">  Developer Note:</span>
<span class="cm">  `idx` cannot be const because may be passed to binary viewer where temporary byte swapping may be done</span>

<span class="cm">.seealso: `<a href="../manualpages/Viewer/PetscViewer.html">PetscViewer</a>`, `<a href="../manualpages/Sys/PetscIntView.html">PetscIntView</a>()`, `<a href="../manualpages/Sys/PetscRealView.html">PetscRealView</a>()`</span>
<span class="cm">@*/</span>
<span class="n"><a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a></span><span class="w"> </span><span class="nf"><a href="../manualpages/Sys/PetscIntViewNumColumns.html">PetscIntViewNumColumns</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscInt.html">PetscInt</a></span><span class="w"> </span><span class="n">N</span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscInt.html">PetscInt</a></span><span class="w"> </span><span class="n">Ncol</span><span class="p">,</span><span class="w"> </span><span class="k">const</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscInt.html">PetscInt</a></span><span class="w"> </span><span class="n">idx</span><span class="p">[],</span><span class="w"> </span><span class="n"><a href="../manualpages/Viewer/PetscViewer.html">PetscViewer</a></span><span class="w"> </span><span class="n">viewer</span><span class="p">)</span>
<span class="p">{</span>
<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscMPIInt.html">PetscMPIInt</a></span><span class="w"> </span><span class="n">rank</span><span class="p">,</span><span class="w"> </span><span class="n">size</span><span class="p">;</span>
<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscInt.html">PetscInt</a></span><span class="w">    </span><span class="n">j</span><span class="p">,</span><span class="w"> </span><span class="n">i</span><span class="p">,</span><span class="w"> </span><span class="n">n</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="n">N</span><span class="w"> </span><span class="o">/</span><span class="w"> </span><span class="n">Ncol</span><span class="p">,</span><span class="w"> </span><span class="n">p</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="n">N</span><span class="w"> </span><span class="o">%</span><span class="w"> </span><span class="n">Ncol</span><span class="p">;</span>
<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscBool.html">PetscBool</a></span><span class="w">   </span><span class="n">iascii</span><span class="p">,</span><span class="w"> </span><span class="n">isbinary</span><span class="p">;</span>
<span class="w">  </span><span class="n"><a href="../manualpages/Sys/MPI_Comm.html">MPI_Comm</a></span><span class="w">    </span><span class="n">comm</span><span class="p">;</span>

<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscFunctionBegin.html">PetscFunctionBegin</a></span><span class="p">;</span>
<span class="w">  </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="o">!</span><span class="n">viewer</span><span class="p">)</span><span class="w"> </span><span class="n">viewer</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="n"><a href="../manualpages/Viewer/PETSC_VIEWER_STDOUT_SELF.html">PETSC_VIEWER_STDOUT_SELF</a></span><span class="p">;</span>
<span class="w">  </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">N</span><span class="p">)</span><span class="w"> </span><span class="n">PetscAssertPointer</span><span class="p">(</span><span class="n">idx</span><span class="p">,</span><span class="w"> </span><span class="mi">3</span><span class="p">);</span>
<span class="w">  </span><span class="n">PetscValidHeaderSpecific</span><span class="p">(</span><span class="n">viewer</span><span class="p">,</span><span class="w"> </span><span class="n">PETSC_VIEWER_CLASSID</span><span class="p">,</span><span class="w"> </span><span class="mi">4</span><span class="p">);</span>
<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscObjectGetComm.html">PetscObjectGetComm</a></span><span class="p">((</span><span class="n"><a href="../manualpages/Sys/PetscObject.html">PetscObject</a></span><span class="p">)</span><span class="n">viewer</span><span class="p">,</span><span class="w"> </span><span class="o">&amp;</span><span class="n">comm</span><span class="p">));</span>
<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscCallMPI.html">PetscCallMPI</a></span><span class="p">(</span><span class="n"><a href="http://www.mpich.org/static/docs/latest/www3/MPI_Comm_size.html#MPI_Comm_size">MPI_Comm_size</a></span><span class="p">(</span><span class="n">comm</span><span class="p">,</span><span class="w"> </span><span class="o">&amp;</span><span class="n">size</span><span class="p">));</span>
<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscCallMPI.html">PetscCallMPI</a></span><span class="p">(</span><span class="n"><a href="http://www.mpich.org/static/docs/latest/www3/MPI_Comm_rank.html#MPI_Comm_rank">MPI_Comm_rank</a></span><span class="p">(</span><span class="n">comm</span><span class="p">,</span><span class="w"> </span><span class="o">&amp;</span><span class="n">rank</span><span class="p">));</span>

<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscObjectTypeCompare.html">PetscObjectTypeCompare</a></span><span class="p">((</span><span class="n"><a href="../manualpages/Sys/PetscObject.html">PetscObject</a></span><span class="p">)</span><span class="n">viewer</span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../manualpages/Viewer/PETSCVIEWERASCII.html">PETSCVIEWERASCII</a></span><span class="p">,</span><span class="w"> </span><span class="o">&amp;</span><span class="n">iascii</span><span class="p">));</span>
<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscObjectTypeCompare.html">PetscObjectTypeCompare</a></span><span class="p">((</span><span class="n"><a href="../manualpages/Sys/PetscObject.html">PetscObject</a></span><span class="p">)</span><span class="n">viewer</span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../manualpages/Viewer/PETSCVIEWERBINARY.html">PETSCVIEWERBINARY</a></span><span class="p">,</span><span class="w"> </span><span class="o">&amp;</span><span class="n">isbinary</span><span class="p">));</span>
<span class="w">  </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">iascii</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">    </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Viewer/PetscViewerASCIIPushSynchronized.html">PetscViewerASCIIPushSynchronized</a></span><span class="p">(</span><span class="n">viewer</span><span class="p">));</span>
<span class="w">    </span><span class="k">for</span><span class="w"> </span><span class="p">(</span><span class="n">i</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="mi">0</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="w"> </span><span class="o">&lt;</span><span class="w"> </span><span class="n">n</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="o">++</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">      </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">size</span><span class="w"> </span><span class="o">&gt;</span><span class="w"> </span><span class="mi">1</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">        </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Viewer/PetscViewerASCIISynchronizedPrintf.html">PetscViewerASCIISynchronizedPrintf</a></span><span class="p">(</span><span class="n">viewer</span><span class="p">,</span><span class="w"> </span><span class="s">&quot;[%d] %&quot;</span><span class="w"> </span><span class="n">PetscInt_FMT</span><span class="w"> </span><span class="s">&quot;:&quot;</span><span class="p">,</span><span class="w"> </span><span class="n">rank</span><span class="p">,</span><span class="w"> </span><span class="n">Ncol</span><span class="w"> </span><span class="o">*</span><span class="w"> </span><span class="n">i</span><span class="p">));</span>
<span class="w">      </span><span class="p">}</span><span class="w"> </span><span class="k">else</span><span class="w"> </span><span class="p">{</span>
<span class="w">        </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Viewer/PetscViewerASCIISynchronizedPrintf.html">PetscViewerASCIISynchronizedPrintf</a></span><span class="p">(</span><span class="n">viewer</span><span class="p">,</span><span class="w"> </span><span class="s">&quot;%&quot;</span><span class="w"> </span><span class="n">PetscInt_FMT</span><span class="w"> </span><span class="s">&quot;:&quot;</span><span class="p">,</span><span class="w"> </span><span class="n">Ncol</span><span class="w"> </span><span class="o">*</span><span class="w"> </span><span class="n">i</span><span class="p">));</span>
<span class="w">      </span><span class="p">}</span>
<span class="w">      </span><span class="k">for</span><span class="w"> </span><span class="p">(</span><span class="n">j</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="mi">0</span><span class="p">;</span><span class="w"> </span><span class="n">j</span><span class="w"> </span><span class="o">&lt;</span><span class="w"> </span><span class="n">Ncol</span><span class="p">;</span><span class="w"> </span><span class="n">j</span><span class="o">++</span><span class="p">)</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Viewer/PetscViewerASCIISynchronizedPrintf.html">PetscViewerASCIISynchronizedPrintf</a></span><span class="p">(</span><span class="n">viewer</span><span class="p">,</span><span class="w"> </span><span class="s">&quot; %&quot;</span><span class="w"> </span><span class="n">PetscInt_FMT</span><span class="p">,</span><span class="w"> </span><span class="n">idx</span><span class="p">[</span><span class="n">i</span><span class="w"> </span><span class="o">*</span><span class="w"> </span><span class="n">Ncol</span><span class="w"> </span><span class="o">+</span><span class="w"> </span><span class="n">j</span><span class="p">]));</span>
<span class="w">      </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Viewer/PetscViewerASCIISynchronizedPrintf.html">PetscViewerASCIISynchronizedPrintf</a></span><span class="p">(</span><span class="n">viewer</span><span class="p">,</span><span class="w"> </span><span class="s">&quot;</span><span class="se">\n</span><span class="s">&quot;</span><span class="p">));</span>
<span class="w">    </span><span class="p">}</span>
<span class="w">    </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">p</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">      </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">size</span><span class="w"> </span><span class="o">&gt;</span><span class="w"> </span><span class="mi">1</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">        </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Viewer/PetscViewerASCIISynchronizedPrintf.html">PetscViewerASCIISynchronizedPrintf</a></span><span class="p">(</span><span class="n">viewer</span><span class="p">,</span><span class="w"> </span><span class="s">&quot;[%d] %&quot;</span><span class="w"> </span><span class="n">PetscInt_FMT</span><span class="w"> </span><span class="s">&quot;:&quot;</span><span class="p">,</span><span class="w"> </span><span class="n">rank</span><span class="p">,</span><span class="w"> </span><span class="n">Ncol</span><span class="w"> </span><span class="o">*</span><span class="w"> </span><span class="n">n</span><span class="p">));</span>
<span class="w">      </span><span class="p">}</span><span class="w"> </span><span class="k">else</span><span class="w"> </span><span class="p">{</span>
<span class="w">        </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Viewer/PetscViewerASCIISynchronizedPrintf.html">PetscViewerASCIISynchronizedPrintf</a></span><span class="p">(</span><span class="n">viewer</span><span class="p">,</span><span class="w"> </span><span class="s">&quot;%&quot;</span><span class="w"> </span><span class="n">PetscInt_FMT</span><span class="w"> </span><span class="s">&quot;:&quot;</span><span class="p">,</span><span class="w"> </span><span class="n">Ncol</span><span class="w"> </span><span class="o">*</span><span class="w"> </span><span class="n">n</span><span class="p">));</span>
<span class="w">      </span><span class="p">}</span>
<span class="w">      </span><span class="k">for</span><span class="w"> </span><span class="p">(</span><span class="n">i</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="mi">0</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="w"> </span><span class="o">&lt;</span><span class="w"> </span><span class="n">p</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="o">++</span><span class="p">)</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Viewer/PetscViewerASCIISynchronizedPrintf.html">PetscViewerASCIISynchronizedPrintf</a></span><span class="p">(</span><span class="n">viewer</span><span class="p">,</span><span class="w"> </span><span class="s">&quot; %&quot;</span><span class="w"> </span><span class="n">PetscInt_FMT</span><span class="p">,</span><span class="w"> </span><span class="n">idx</span><span class="p">[</span><span class="n">Ncol</span><span class="w"> </span><span class="o">*</span><span class="w"> </span><span class="n">n</span><span class="w"> </span><span class="o">+</span><span class="w"> </span><span class="n">i</span><span class="p">]));</span>
<span class="w">      </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Viewer/PetscViewerASCIISynchronizedPrintf.html">PetscViewerASCIISynchronizedPrintf</a></span><span class="p">(</span><span class="n">viewer</span><span class="p">,</span><span class="w"> </span><span class="s">&quot;</span><span class="se">\n</span><span class="s">&quot;</span><span class="p">));</span>
<span class="w">    </span><span class="p">}</span>
<span class="w">    </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Viewer/PetscViewerFlush.html">PetscViewerFlush</a></span><span class="p">(</span><span class="n">viewer</span><span class="p">));</span>
<span class="w">    </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Viewer/PetscViewerASCIIPopSynchronized.html">PetscViewerASCIIPopSynchronized</a></span><span class="p">(</span><span class="n">viewer</span><span class="p">));</span>
<span class="w">  </span><span class="p">}</span><span class="w"> </span><span class="k">else</span><span class="w"> </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">isbinary</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">    </span><span class="n"><a href="../manualpages/Sys/PetscMPIInt.html">PetscMPIInt</a></span><span class="w"> </span><span class="o">*</span><span class="n">sizes</span><span class="p">,</span><span class="w"> </span><span class="n">Ntotal</span><span class="p">,</span><span class="w"> </span><span class="o">*</span><span class="n">displs</span><span class="p">,</span><span class="w"> </span><span class="n">NN</span><span class="p">;</span>
<span class="w">    </span><span class="n"><a href="../manualpages/Sys/PetscInt.html">PetscInt</a></span><span class="w">    </span><span class="o">*</span><span class="n">array</span><span class="p">;</span>

<span class="w">    </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscMPIIntCast.html">PetscMPIIntCast</a></span><span class="p">(</span><span class="n">N</span><span class="p">,</span><span class="w"> </span><span class="o">&amp;</span><span class="n">NN</span><span class="p">));</span>

<span class="w">    </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">size</span><span class="w"> </span><span class="o">&gt;</span><span class="w"> </span><span class="mi">1</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">      </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">rank</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">        </span><span class="n"><a href="../manualpages/Sys/PetscCallMPI.html">PetscCallMPI</a></span><span class="p">(</span><span class="n"><a href="http://www.mpich.org/static/docs/latest/www3/MPI_Gather.html#MPI_Gather">MPI_Gather</a></span><span class="p">(</span><span class="o">&amp;</span><span class="n">NN</span><span class="p">,</span><span class="w"> </span><span class="mi">1</span><span class="p">,</span><span class="w"> </span><span class="n">MPI_INT</span><span class="p">,</span><span class="w"> </span><span class="nb">NULL</span><span class="p">,</span><span class="w"> </span><span class="mi">0</span><span class="p">,</span><span class="w"> </span><span class="n">MPI_INT</span><span class="p">,</span><span class="w"> </span><span class="mi">0</span><span class="p">,</span><span class="w"> </span><span class="n">comm</span><span class="p">));</span>
<span class="w">        </span><span class="n"><a href="../manualpages/Sys/PetscCallMPI.html">PetscCallMPI</a></span><span class="p">(</span><span class="n"><a href="http://www.mpich.org/static/docs/latest/www3/MPI_Gatherv.html#MPI_Gatherv">MPI_Gatherv</a></span><span class="p">(</span><span class="n">idx</span><span class="p">,</span><span class="w"> </span><span class="n">NN</span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/MPIU_INT.html">MPIU_INT</a></span><span class="p">,</span><span class="w"> </span><span class="nb">NULL</span><span class="p">,</span><span class="w"> </span><span class="nb">NULL</span><span class="p">,</span><span class="w"> </span><span class="nb">NULL</span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/MPIU_INT.html">MPIU_INT</a></span><span class="p">,</span><span class="w"> </span><span class="mi">0</span><span class="p">,</span><span class="w"> </span><span class="n">comm</span><span class="p">));</span>
<span class="w">      </span><span class="p">}</span><span class="w"> </span><span class="k">else</span><span class="w"> </span><span class="p">{</span>
<span class="w">        </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscMalloc1.html">PetscMalloc1</a></span><span class="p">(</span><span class="n">size</span><span class="p">,</span><span class="w"> </span><span class="o">&amp;</span><span class="n">sizes</span><span class="p">));</span>
<span class="w">        </span><span class="n"><a href="../manualpages/Sys/PetscCallMPI.html">PetscCallMPI</a></span><span class="p">(</span><span class="n"><a href="http://www.mpich.org/static/docs/latest/www3/MPI_Gather.html#MPI_Gather">MPI_Gather</a></span><span class="p">(</span><span class="o">&amp;</span><span class="n">NN</span><span class="p">,</span><span class="w"> </span><span class="mi">1</span><span class="p">,</span><span class="w"> </span><span class="n">MPI_INT</span><span class="p">,</span><span class="w"> </span><span class="n">sizes</span><span class="p">,</span><span class="w"> </span><span class="mi">1</span><span class="p">,</span><span class="w"> </span><span class="n">MPI_INT</span><span class="p">,</span><span class="w"> </span><span class="mi">0</span><span class="p">,</span><span class="w"> </span><span class="n">comm</span><span class="p">));</span>
<span class="w">        </span><span class="n">Ntotal</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="n">sizes</span><span class="p">[</span><span class="mi">0</span><span class="p">];</span>
<span class="w">        </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscMalloc1.html">PetscMalloc1</a></span><span class="p">(</span><span class="n">size</span><span class="p">,</span><span class="w"> </span><span class="o">&amp;</span><span class="n">displs</span><span class="p">));</span>
<span class="w">        </span><span class="n">displs</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="mi">0</span><span class="p">;</span>
<span class="w">        </span><span class="k">for</span><span class="w"> </span><span class="p">(</span><span class="n">i</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="mi">1</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="w"> </span><span class="o">&lt;</span><span class="w"> </span><span class="n">size</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="o">++</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">          </span><span class="n">Ntotal</span><span class="w"> </span><span class="o">+=</span><span class="w"> </span><span class="n">sizes</span><span class="p">[</span><span class="n">i</span><span class="p">];</span>
<span class="w">          </span><span class="n">displs</span><span class="p">[</span><span class="n">i</span><span class="p">]</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="n">displs</span><span class="p">[</span><span class="n">i</span><span class="w"> </span><span class="o">-</span><span class="w"> </span><span class="mi">1</span><span class="p">]</span><span class="w"> </span><span class="o">+</span><span class="w"> </span><span class="n">sizes</span><span class="p">[</span><span class="n">i</span><span class="w"> </span><span class="o">-</span><span class="w"> </span><span class="mi">1</span><span class="p">];</span>
<span class="w">        </span><span class="p">}</span>
<span class="w">        </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscMalloc1.html">PetscMalloc1</a></span><span class="p">(</span><span class="n">Ntotal</span><span class="p">,</span><span class="w"> </span><span class="o">&amp;</span><span class="n">array</span><span class="p">));</span>
<span class="w">        </span><span class="n"><a href="../manualpages/Sys/PetscCallMPI.html">PetscCallMPI</a></span><span class="p">(</span><span class="n"><a href="http://www.mpich.org/static/docs/latest/www3/MPI_Gatherv.html#MPI_Gatherv">MPI_Gatherv</a></span><span class="p">(</span><span class="n">idx</span><span class="p">,</span><span class="w"> </span><span class="n">NN</span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/MPIU_INT.html">MPIU_INT</a></span><span class="p">,</span><span class="w"> </span><span class="n">array</span><span class="p">,</span><span class="w"> </span><span class="n">sizes</span><span class="p">,</span><span class="w"> </span><span class="n">displs</span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/MPIU_INT.html">MPIU_INT</a></span><span class="p">,</span><span class="w"> </span><span class="mi">0</span><span class="p">,</span><span class="w"> </span><span class="n">comm</span><span class="p">));</span>
<span class="w">        </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Viewer/PetscViewerBinaryWrite.html">PetscViewerBinaryWrite</a></span><span class="p">(</span><span class="n">viewer</span><span class="p">,</span><span class="w"> </span><span class="n">array</span><span class="p">,</span><span class="w"> </span><span class="n">Ntotal</span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscDataType.html">PETSC_INT</a></span><span class="p">));</span>
<span class="w">        </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscFree.html">PetscFree</a></span><span class="p">(</span><span class="n">sizes</span><span class="p">));</span>
<span class="w">        </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscFree.html">PetscFree</a></span><span class="p">(</span><span class="n">displs</span><span class="p">));</span>
<span class="w">        </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscFree.html">PetscFree</a></span><span class="p">(</span><span class="n">array</span><span class="p">));</span>
<span class="w">      </span><span class="p">}</span>
<span class="w">    </span><span class="p">}</span><span class="w"> </span><span class="k">else</span><span class="w"> </span><span class="p">{</span>
<span class="w">      </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Viewer/PetscViewerBinaryWrite.html">PetscViewerBinaryWrite</a></span><span class="p">(</span><span class="n">viewer</span><span class="p">,</span><span class="w"> </span><span class="n">idx</span><span class="p">,</span><span class="w"> </span><span class="n">N</span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscDataType.html">PETSC_INT</a></span><span class="p">));</span>
<span class="w">    </span><span class="p">}</span>
<span class="w">  </span><span class="p">}</span><span class="w"> </span><span class="k">else</span><span class="w"> </span><span class="p">{</span>
<span class="w">    </span><span class="k">const</span><span class="w"> </span><span class="kt">char</span><span class="w"> </span><span class="o">*</span><span class="n">tname</span><span class="p">;</span>
<span class="w">    </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscObjectGetName.html">PetscObjectGetName</a></span><span class="p">((</span><span class="n"><a href="../manualpages/Sys/PetscObject.html">PetscObject</a></span><span class="p">)</span><span class="n">viewer</span><span class="p">,</span><span class="w"> </span><span class="o">&amp;</span><span class="n">tname</span><span class="p">));</span>
<span class="w">    </span><span class="n"><a href="../manualpages/Sys/SETERRQ.html">SETERRQ</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PETSC_COMM_SELF.html">PETSC_COMM_SELF</a></span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscErrorCode.html">PETSC_ERR_SUP</a></span><span class="p">,</span><span class="w"> </span><span class="s">&quot;Cannot handle that <a href="../manualpages/Viewer/PetscViewer.html">PetscViewer</a> of type %s&quot;</span><span class="p">,</span><span class="w"> </span><span class="n">tname</span><span class="p">);</span>
<span class="w">  </span><span class="p">}</span>
<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscFunctionReturn.html">PetscFunctionReturn</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscErrorCode.html">PETSC_SUCCESS</a></span><span class="p">);</span>
<span class="p">}</span>
</pre></div>
</div>
<p>You can set a new error handler with the command
<code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscPushErrorHandler.html">PetscPushErrorHandler</a>()</span></code>, which maintains a linked list of error
handlers. The most recent error handler is removed via
<code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscPopErrorHandler.html">PetscPopErrorHandler</a>()</span></code>.</p>
<p>PETSc provides several default error handlers:</p>
<ul class="simple">
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscTraceBackErrorHandler.html">PetscTraceBackErrorHandler</a>()</span></code>, the default;</p></li>
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscAbortErrorHandler.html">PetscAbortErrorHandler</a>()</span></code>, called with <code class="docutils notranslate"><span class="pre">-onerrorabort</span></code>, useful when running in the debugger;</p></li>
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscReturnErrorHandler.html">PetscReturnErrorHandler</a>()</span></code>, which returns up the stack without printing error messages;</p></li>
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscEmacsClientErrorHandler.html">PetscEmacsClientErrorHandler</a>()</span></code>;</p></li>
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscMPIAbortErrorHandler.html">PetscMPIAbortErrorHandler</a>()</span></code>, which calls <code class="docutils notranslate"><span class="pre">MPIAbort()</span></code> after printing the error message; and</p></li>
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscAttachDebuggerErrorHandler.html">PetscAttachDebuggerErrorHandler</a>()</span></code>, called with <code class="docutils notranslate"><span class="pre">-onerrorattachdebugger</span></code>.</p></li>
</ul>
</section>
<section id="error-codes">
<h3>Error Codes<a class="headerlink" href="#error-codes" title="Link to this heading">#</a></h3>
<p>The PETSc error handler takes an error code. The generic error codes are
defined in
<a href="../include/petscerror.h.html">include/petscerror.h</a>
The same error code is used many times in the libraries. For example,
the error code <code class="docutils notranslate"><span class="pre">PETSCERRMEM</span></code> is used whenever a requested memory
allocation is not available.</p>
</section>
<section id="detailed-error-messages">
<h3>Detailed Error Messages<a class="headerlink" href="#detailed-error-messages" title="Link to this heading">#</a></h3>
<p>In a modern parallel component-oriented application code, it does not
always make sense to simply print error messages to the terminal (and
more than likely there is no “terminal”, for example, with Microsoft
Windows or Apple iPad applications). PETSc provides the replaceable
function pointer</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="p">(</span><span class="o">*</span><span class="n"><a href="../manualpages/Sys/PetscErrorPrintf.html">PetscErrorPrintf</a></span><span class="p">)(</span><span class="s">&quot;Format&quot;</span><span class="p">,...);</span>
</pre></div>
</div>
<p>which, by default, prints to standard out. Thus, error messages should
not be printed with <code class="docutils notranslate"><span class="pre">printf()</span></code> or <code class="docutils notranslate"><span class="pre">fprintf()</span></code>. Rather, they should
be printed with <code class="docutils notranslate"><span class="pre">(*<a href="../manualpages/Sys/PetscErrorPrintf.html">PetscErrorPrintf</a>)()</span></code>. You can direct all error
messages to <code class="docutils notranslate"><span class="pre">stderr</span></code>, instead of the default <code class="docutils notranslate"><span class="pre">stdout</span></code>, with the
command line option <code class="docutils notranslate"><span class="pre">-erroroutputstderr</span></code>.</p>
</section>
<section id="c-exceptions">
<h3>C++ Exceptions<a class="headerlink" href="#c-exceptions" title="Link to this heading">#</a></h3>
<p>In PETSc code, when one calls C++ functions that do not return with an error code but might
instead throw C++ exceptions, one can use <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/CHKERRCXX.html">CHKERRCXX</a>(func)</span></code>, which catches the exceptions
in <em>func</em> and then calls <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/SETERRQ.html">SETERRQ</a>()</span></code>. The macro <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/CHKERRCXX.html">CHKERRCXX</a>(func)</span></code> is given by</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="cp">#define <a href="../manualpages/Sys/CHKERRCXX.html">CHKERRCXX</a>(...) <a href="../manualpages/Sys/PetscCallCXX.html">PetscCallCXX</a>(__VA_ARGS__)</span>
</pre></div>
</div>
</section>
</section>
<section id="memory-management">
<h2>Memory Management<a class="headerlink" href="#memory-management" title="Link to this heading">#</a></h2>
<p>PETSc provides simple wrappers for the system <code class="docutils notranslate"><span class="pre">malloc(),</span> <span class="pre">calloc()</span></code>,
and <code class="docutils notranslate"><span class="pre">free()</span></code> routines. The public interface for these is provided in
<code class="docutils notranslate"><span class="pre">petscsys.h</span></code>, while the implementation code is in <code class="docutils notranslate"><span class="pre">src/sys/memory</span></code>.
The most basic interfaces are</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="cp">#define <a href="../manualpages/Sys/PetscMalloc.html">PetscMalloc</a>(a, b) ((*PetscTrMalloc)((a), <a href="../manualpages/Sys/PETSC_FALSE.html">PETSC_FALSE</a>, __LINE__, PETSC_FUNCTION_NAME, __FILE__, (void **)(b)))</span>
</pre></div>
</div>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="cp">#define <a href="../manualpages/Sys/PetscFree.html">PetscFree</a>(a) ((<a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a>)((*PetscTrFree)((void *)(a), __LINE__, PETSC_FUNCTION_NAME, __FILE__) || ((a) = <a href="../manualpages/Sys/PETSC_NULLPTR.html">PETSC_NULLPTR</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PETSC_SUCCESS</a>)))</span>
</pre></div>
</div>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a></span><span class="w"> </span><span class="nf"><a href="../manualpages/Sys/PetscMallocA.html">PetscMallocA</a></span><span class="p">(</span><span class="kt">int</span><span class="w"> </span><span class="n">n</span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscBool.html">PetscBool</a></span><span class="w"> </span><span class="n">clear</span><span class="p">,</span><span class="w"> </span><span class="kt">int</span><span class="w"> </span><span class="n">lineno</span><span class="p">,</span><span class="w"> </span><span class="k">const</span><span class="w"> </span><span class="kt">char</span><span class="w"> </span><span class="o">*</span><span class="n">function</span><span class="p">,</span><span class="w"> </span><span class="k">const</span><span class="w"> </span><span class="kt">char</span><span class="w"> </span><span class="o">*</span><span class="n">filename</span><span class="p">,</span><span class="w"> </span><span class="kt">size_t</span><span class="w"> </span><span class="n">bytes0</span><span class="p">,</span><span class="w"> </span><span class="kt">void</span><span class="w"> </span><span class="o">*</span><span class="n">ptr0</span><span class="p">,</span><span class="w"> </span><span class="p">...)</span>
<span class="p">{</span>
<span class="w">  </span><span class="kt">va_list</span><span class="w"> </span><span class="n">Argp</span><span class="p">;</span>
<span class="w">  </span><span class="kt">size_t</span><span class="w">  </span><span class="n">bytes</span><span class="p">[</span><span class="mi">8</span><span class="p">],</span><span class="w"> </span><span class="n">sumbytes</span><span class="p">;</span>
<span class="w">  </span><span class="kt">void</span><span class="w">  </span><span class="o">**</span><span class="n">ptr</span><span class="p">[</span><span class="mi">8</span><span class="p">];</span>
<span class="w">  </span><span class="kt">int</span><span class="w">     </span><span class="n">i</span><span class="p">;</span>

<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscFunctionBegin.html">PetscFunctionBegin</a></span><span class="p">;</span>
<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscCheck.html">PetscCheck</a></span><span class="p">(</span><span class="n">n</span><span class="w"> </span><span class="o">&lt;=</span><span class="w"> </span><span class="mi">8</span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PETSC_COMM_SELF.html">PETSC_COMM_SELF</a></span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscErrorCode.html">PETSC_ERR_ARG_OUTOFRANGE</a></span><span class="p">,</span><span class="w"> </span><span class="s">&quot;Attempt to allocate %d objects but only 8 supported&quot;</span><span class="p">,</span><span class="w"> </span><span class="n">n</span><span class="p">);</span>
<span class="w">  </span><span class="n">bytes</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="n">bytes0</span><span class="p">;</span>
<span class="w">  </span><span class="n">ptr</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span><span class="w">   </span><span class="o">=</span><span class="w"> </span><span class="p">(</span><span class="kt">void</span><span class="w"> </span><span class="o">**</span><span class="p">)</span><span class="n">ptr0</span><span class="p">;</span>
<span class="w">  </span><span class="n">sumbytes</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="p">(</span><span class="n">bytes0</span><span class="w"> </span><span class="o">+</span><span class="w"> </span><span class="n">PETSC_MEMALIGN</span><span class="w"> </span><span class="o">-</span><span class="w"> </span><span class="mi">1</span><span class="p">)</span><span class="w"> </span><span class="o">&amp;</span><span class="w"> </span><span class="o">~</span><span class="p">(</span><span class="n">PETSC_MEMALIGN</span><span class="w"> </span><span class="o">-</span><span class="w"> </span><span class="mi">1</span><span class="p">);</span>
<span class="w">  </span><span class="n">va_start</span><span class="p">(</span><span class="n">Argp</span><span class="p">,</span><span class="w"> </span><span class="n">ptr0</span><span class="p">);</span>
<span class="w">  </span><span class="k">for</span><span class="w"> </span><span class="p">(</span><span class="n">i</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="mi">1</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="w"> </span><span class="o">&lt;</span><span class="w"> </span><span class="n">n</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="o">++</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">    </span><span class="n">bytes</span><span class="p">[</span><span class="n">i</span><span class="p">]</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="n">va_arg</span><span class="p">(</span><span class="n">Argp</span><span class="p">,</span><span class="w"> </span><span class="kt">size_t</span><span class="p">);</span>
<span class="w">    </span><span class="n">ptr</span><span class="p">[</span><span class="n">i</span><span class="p">]</span><span class="w">   </span><span class="o">=</span><span class="w"> </span><span class="n">va_arg</span><span class="p">(</span><span class="n">Argp</span><span class="p">,</span><span class="w"> </span><span class="kt">void</span><span class="w"> </span><span class="o">**</span><span class="p">);</span>
<span class="w">    </span><span class="n">sumbytes</span><span class="w"> </span><span class="o">+=</span><span class="w"> </span><span class="p">(</span><span class="n">bytes</span><span class="p">[</span><span class="n">i</span><span class="p">]</span><span class="w"> </span><span class="o">+</span><span class="w"> </span><span class="n">PETSC_MEMALIGN</span><span class="w"> </span><span class="o">-</span><span class="w"> </span><span class="mi">1</span><span class="p">)</span><span class="w"> </span><span class="o">&amp;</span><span class="w"> </span><span class="o">~</span><span class="p">(</span><span class="n">PETSC_MEMALIGN</span><span class="w"> </span><span class="o">-</span><span class="w"> </span><span class="mi">1</span><span class="p">);</span>
<span class="w">  </span><span class="p">}</span>
<span class="w">  </span><span class="n">va_end</span><span class="p">(</span><span class="n">Argp</span><span class="p">);</span>
<span class="w">  </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">petscmalloccoalesce</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">    </span><span class="kt">char</span><span class="w"> </span><span class="o">*</span><span class="n">p</span><span class="p">;</span>
<span class="w">    </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">((</span><span class="o">*</span><span class="n">PetscTrMalloc</span><span class="p">)(</span><span class="n">sumbytes</span><span class="p">,</span><span class="w"> </span><span class="n">clear</span><span class="p">,</span><span class="w"> </span><span class="n">lineno</span><span class="p">,</span><span class="w"> </span><span class="n">function</span><span class="p">,</span><span class="w"> </span><span class="n">filename</span><span class="p">,</span><span class="w"> </span><span class="p">(</span><span class="kt">void</span><span class="w"> </span><span class="o">**</span><span class="p">)</span><span class="o">&amp;</span><span class="n">p</span><span class="p">));</span>
<span class="w">    </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">p</span><span class="w"> </span><span class="o">==</span><span class="w"> </span><span class="nb">NULL</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">      </span><span class="k">for</span><span class="w"> </span><span class="p">(</span><span class="n">i</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="mi">0</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="w"> </span><span class="o">&lt;</span><span class="w"> </span><span class="n">n</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="o">++</span><span class="p">)</span><span class="w"> </span><span class="o">*</span><span class="n">ptr</span><span class="p">[</span><span class="n">i</span><span class="p">]</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="nb">NULL</span><span class="p">;</span>
<span class="w">    </span><span class="p">}</span><span class="w"> </span><span class="k">else</span><span class="w"> </span><span class="p">{</span>
<span class="w">      </span><span class="k">for</span><span class="w"> </span><span class="p">(</span><span class="n">i</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="mi">0</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="w"> </span><span class="o">&lt;</span><span class="w"> </span><span class="n">n</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="o">++</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">        </span><span class="o">*</span><span class="n">ptr</span><span class="p">[</span><span class="n">i</span><span class="p">]</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="n">bytes</span><span class="p">[</span><span class="n">i</span><span class="p">]</span><span class="w"> </span><span class="o">?</span><span class="w"> </span><span class="n">p</span><span class="w"> </span><span class="o">:</span><span class="w"> </span><span class="nb">NULL</span><span class="p">;</span>
<span class="w">        </span><span class="n">p</span><span class="w">       </span><span class="o">=</span><span class="w"> </span><span class="p">(</span><span class="kt">char</span><span class="w"> </span><span class="o">*</span><span class="p">)</span><span class="n"><a href="../manualpages/Sys/PetscAddrAlign.html">PetscAddrAlign</a></span><span class="p">(</span><span class="n">p</span><span class="w"> </span><span class="o">+</span><span class="w"> </span><span class="n">bytes</span><span class="p">[</span><span class="n">i</span><span class="p">]);</span>
<span class="w">      </span><span class="p">}</span>
<span class="w">    </span><span class="p">}</span>
<span class="w">  </span><span class="p">}</span><span class="w"> </span><span class="k">else</span><span class="w"> </span><span class="p">{</span>
<span class="w">    </span><span class="k">for</span><span class="w"> </span><span class="p">(</span><span class="n">i</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="mi">0</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="w"> </span><span class="o">&lt;</span><span class="w"> </span><span class="n">n</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="o">++</span><span class="p">)</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">((</span><span class="o">*</span><span class="n">PetscTrMalloc</span><span class="p">)(</span><span class="n">bytes</span><span class="p">[</span><span class="n">i</span><span class="p">],</span><span class="w"> </span><span class="n">clear</span><span class="p">,</span><span class="w"> </span><span class="n">lineno</span><span class="p">,</span><span class="w"> </span><span class="n">function</span><span class="p">,</span><span class="w"> </span><span class="n">filename</span><span class="p">,</span><span class="w"> </span><span class="n">ptr</span><span class="p">[</span><span class="n">i</span><span class="p">]));</span>
<span class="w">  </span><span class="p">}</span>
<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscFunctionReturn.html">PetscFunctionReturn</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscErrorCode.html">PETSC_SUCCESS</a></span><span class="p">);</span>
<span class="p">}</span>
</pre></div>
</div>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a></span><span class="w"> </span><span class="nf"><a href="../manualpages/Sys/PetscFreeA.html">PetscFreeA</a></span><span class="p">(</span><span class="kt">int</span><span class="w"> </span><span class="n">n</span><span class="p">,</span><span class="w"> </span><span class="kt">int</span><span class="w"> </span><span class="n">lineno</span><span class="p">,</span><span class="w"> </span><span class="k">const</span><span class="w"> </span><span class="kt">char</span><span class="w"> </span><span class="o">*</span><span class="n">function</span><span class="p">,</span><span class="w"> </span><span class="k">const</span><span class="w"> </span><span class="kt">char</span><span class="w"> </span><span class="o">*</span><span class="n">filename</span><span class="p">,</span><span class="w"> </span><span class="kt">void</span><span class="w"> </span><span class="o">*</span><span class="n">ptr0</span><span class="p">,</span><span class="w"> </span><span class="p">...)</span>
<span class="p">{</span>
<span class="w">  </span><span class="kt">va_list</span><span class="w"> </span><span class="n">Argp</span><span class="p">;</span>
<span class="w">  </span><span class="kt">void</span><span class="w">  </span><span class="o">**</span><span class="n">ptr</span><span class="p">[</span><span class="mi">8</span><span class="p">];</span>
<span class="w">  </span><span class="kt">int</span><span class="w">     </span><span class="n">i</span><span class="p">;</span>

<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscFunctionBegin.html">PetscFunctionBegin</a></span><span class="p">;</span>
<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscCheck.html">PetscCheck</a></span><span class="p">((</span><span class="n">n</span><span class="w"> </span><span class="o">&gt;=</span><span class="w"> </span><span class="mi">1</span><span class="p">)</span><span class="w"> </span><span class="o">&amp;&amp;</span><span class="w"> </span><span class="p">(</span><span class="n">n</span><span class="w"> </span><span class="o">&lt;=</span><span class="w"> </span><span class="mi">8</span><span class="p">),</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PETSC_COMM_SELF.html">PETSC_COMM_SELF</a></span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../manualpages/Sys/PetscErrorCode.html">PETSC_ERR_ARG_OUTOFRANGE</a></span><span class="p">,</span><span class="w"> </span><span class="s">&quot;Attempt to allocate %d objects but only up to 8 supported&quot;</span><span class="p">,</span><span class="w"> </span><span class="n">n</span><span class="p">);</span>
<span class="w">  </span><span class="n">ptr</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="p">(</span><span class="kt">void</span><span class="w"> </span><span class="o">**</span><span class="p">)</span><span class="n">ptr0</span><span class="p">;</span>
<span class="w">  </span><span class="n">va_start</span><span class="p">(</span><span class="n">Argp</span><span class="p">,</span><span class="w"> </span><span class="n">ptr0</span><span class="p">);</span>
<span class="w">  </span><span class="k">for</span><span class="w"> </span><span class="p">(</span><span class="n">i</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="mi">1</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="w"> </span><span class="o">&lt;</span><span class="w"> </span><span class="n">n</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="o">++</span><span class="p">)</span><span class="w"> </span><span class="n">ptr</span><span class="p">[</span><span class="n">i</span><span class="p">]</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="n">va_arg</span><span class="p">(</span><span class="n">Argp</span><span class="p">,</span><span class="w"> </span><span class="kt">void</span><span class="w"> </span><span class="o">**</span><span class="p">);</span>
<span class="w">  </span><span class="n">va_end</span><span class="p">(</span><span class="n">Argp</span><span class="p">);</span>
<span class="w">  </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="n">petscmalloccoalesce</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">    </span><span class="k">for</span><span class="w"> </span><span class="p">(</span><span class="n">i</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="mi">0</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="w"> </span><span class="o">&lt;</span><span class="w"> </span><span class="n">n</span><span class="p">;</span><span class="w"> </span><span class="n">i</span><span class="o">++</span><span class="p">)</span><span class="w"> </span><span class="p">{</span><span class="w"> </span><span class="cm">/* Find first nonempty allocation */</span>
<span class="w">      </span><span class="k">if</span><span class="w"> </span><span class="p">(</span><span class="o">*</span><span class="n">ptr</span><span class="p">[</span><span class="n">i</span><span class="p">])</span><span class="w"> </span><span class="k">break</span><span class="p">;</span>
<span class="w">    </span><span class="p">}</span>
<span class="w">    </span><span class="k">while</span><span class="w"> </span><span class="p">(</span><span class="o">--</span><span class="n">n</span><span class="w"> </span><span class="o">&gt;</span><span class="w"> </span><span class="n">i</span><span class="p">)</span><span class="w"> </span><span class="o">*</span><span class="n">ptr</span><span class="p">[</span><span class="n">n</span><span class="p">]</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="nb">NULL</span><span class="p">;</span>
<span class="w">    </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">((</span><span class="o">*</span><span class="n">PetscTrFree</span><span class="p">)(</span><span class="o">*</span><span class="n">ptr</span><span class="p">[</span><span class="n">n</span><span class="p">],</span><span class="w"> </span><span class="n">lineno</span><span class="p">,</span><span class="w"> </span><span class="n">function</span><span class="p">,</span><span class="w"> </span><span class="n">filename</span><span class="p">));</span>
<span class="w">    </span><span class="o">*</span><span class="n">ptr</span><span class="p">[</span><span class="n">n</span><span class="p">]</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="nb">NULL</span><span class="p">;</span>
<span class="w">  </span><span class="p">}</span><span class="w"> </span><span class="k">else</span><span class="w"> </span><span class="p">{</span>
<span class="w">    </span><span class="k">while</span><span class="w"> </span><span class="p">(</span><span class="o">--</span><span class="n">n</span><span class="w"> </span><span class="o">&gt;=</span><span class="w"> </span><span class="mi">0</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="w">      </span><span class="n"><a href="../manualpages/Sys/PetscCall.html">PetscCall</a></span><span class="p">((</span><span class="o">*</span><span class="n">PetscTrFree</span><span class="p">)(</span><span class="o">*</span><span class="n">ptr</span><span class="p">[</span><span class="n">n</span><span class="p">],</span><span class="w"> </span><span class="n">lineno</span><span class="p">,</span><span class="w"> </span><span class="n">function</span><span class="p">,</span><span class="w"> </span><span class="n">filename</span><span class="p">));</span>
<span class="w">      </span><span class="o">*</span><span class="n">ptr</span><span class="p">[</span><span class="n">n</span><span class="p">]</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="nb">NULL</span><span class="p">;</span>
<span class="w">    </span><span class="p">}</span>
<span class="w">  </span><span class="p">}</span>
<span class="w">  </span><span class="n"><a href="../manualpages/Sys/PetscFunctionReturn.html">PetscFunctionReturn</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscErrorCode.html">PETSC_SUCCESS</a></span><span class="p">);</span>
<span class="p">}</span>
</pre></div>
</div>
<p>which allow the use of any number of profiling and error-checking
wrappers for <code class="docutils notranslate"><span class="pre">malloc(),</span> <span class="pre">calloc()</span></code>, and <code class="docutils notranslate"><span class="pre">free()</span></code>. Both
<code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscMallocA.html">PetscMallocA</a>()</span></code> and <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscFreeA.html">PetscFreeA</a>()</span></code> call the function pointer values
<code class="docutils notranslate"><span class="pre">(*PetscTrMalloc)</span></code> and <code class="docutils notranslate"><span class="pre">(*PetscTrFree)</span></code>. <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscMallocSet.html">PetscMallocSet</a>()</span></code> is
used to set these function pointers. The functions are guaranteed to
support requests for zero bytes of memory correctly. Freeing memory
locations also sets the pointer value to zero, preventing later code
from accidentally using memory that has been freed. All PETSc memory
allocation calls are memory aligned on at least double-precision
boundaries; the macro generated by configure <code class="docutils notranslate"><span class="pre">PETSCMEMALIGN</span></code> indicates
in bytes what alignment all allocations have. This can be controlled at
configure time with the option <code class="docutils notranslate"><span class="pre">-with-memalign=&lt;4,8,16,32,64&gt;</span></code>.</p>
<p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscMallocA.html">PetscMallocA</a>()</span></code> supports a request for up to 7 distinct memory
locations of possibly different types. This serves two purposes: it
reduces the number of system <code class="docutils notranslate"><span class="pre">malloc()</span></code> calls, thus potentially
increasing performance, and it clarifies in the code related memory
allocations that should be freed together.</p>
<p>The following macros are the preferred way to obtain and release memory
in the PETSc source code. They automatically manage calling
<code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscMallocA.html">PetscMallocA</a>()</span></code> and <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscFreeA.html">PetscFreeA</a>()</span></code> with the appropriate location
information.</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="cp">#define <a href="../manualpages/Sys/PetscMalloc1.html">PetscMalloc1</a>(m1, r1) <a href="../manualpages/Sys/PetscMallocA.html">PetscMallocA</a>(1, <a href="../manualpages/Sys/PETSC_FALSE.html">PETSC_FALSE</a>, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1))</span>
</pre></div>
</div>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="cp">#define <a href="../manualpages/Sys/PetscMalloc2.html">PetscMalloc2</a>(m1, r1, m2, r2) <a href="../manualpages/Sys/PetscMallocA.html">PetscMallocA</a>(2, <a href="../manualpages/Sys/PETSC_FALSE.html">PETSC_FALSE</a>, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1), ((size_t)((size_t)m2) * sizeof(**(r2))), (r2))</span>
</pre></div>
</div>
<p>…</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="cp">#define <a href="../manualpages/Sys/PetscMalloc7.html">PetscMalloc7</a>(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6, m7, r7) \</span>
</pre></div>
</div>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="cp">#define <a href="../manualpages/Sys/PetscFree.html">PetscFree</a>(a) ((<a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a>)((*PetscTrFree)((void *)(a), __LINE__, PETSC_FUNCTION_NAME, __FILE__) || ((a) = <a href="../manualpages/Sys/PETSC_NULLPTR.html">PETSC_NULLPTR</a>, <a href="../manualpages/Sys/PetscErrorCode.html">PETSC_SUCCESS</a>)))</span>
</pre></div>
</div>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="cp">#define <a href="../manualpages/Sys/PetscFree2.html">PetscFree2</a>(m1, m2) <a href="../manualpages/Sys/PetscFreeA.html">PetscFreeA</a>(2, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &amp;(m1), &amp;(m2))</span>
</pre></div>
</div>
<p>…</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="cp">#define <a href="../manualpages/Sys/PetscFree7.html">PetscFree7</a>(m1, m2, m3, m4, m5, m6, m7) <a href="../manualpages/Sys/PetscFreeA.html">PetscFreeA</a>(7, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &amp;(m1), &amp;(m2), &amp;(m3), &amp;(m4), &amp;(m5), &amp;(m6), &amp;(m7))</span>
</pre></div>
</div>
<p>Similar routines, <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscCalloc1.html">PetscCalloc1</a>()</span></code> to <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscCalloc7.html">PetscCalloc7</a>()</span></code>, provide
memory initialized to zero. The size requests for these macros are in
number of data items requested, not in bytes. This decreases the number
of errors in the code since the compiler determines their sizes from the
object type instead of requiring the user to provide the correct value
with <code class="docutils notranslate"><span class="pre">sizeof()</span></code>.</p>
<p>The routines <code class="docutils notranslate"><span class="pre">PetscTrMallocDefault()</span></code> and <code class="docutils notranslate"><span class="pre">PetscTrFreeDefault()</span></code>,
which are set with the routine <code class="docutils notranslate"><span class="pre">PetscSetUseTrMallocPrivate()</span></code> (and are
used by default for the debug version of PETSc), provide simple logging
and error checking versions of memory allocation.</p>
</section>
<section id="implementation-of-profiling">
<h2>Implementation of Profiling<a class="headerlink" href="#implementation-of-profiling" title="Link to this heading">#</a></h2>
<p>This section provides details about the implementation of event logging
and profiling within the PETSc kernel. The interface for profiling in
PETSc is contained in the file
<a href="../include/petsclog.h.html">include/petsclog.h</a>
The source code for the profile logging is in <code class="docutils notranslate"><span class="pre">src/sys/plog/</span></code>.</p>
<section id="profiling-object-creation-and-destruction">
<h3>Profiling Object Creation and Destruction<a class="headerlink" href="#profiling-object-creation-and-destruction" title="Link to this heading">#</a></h3>
<p>The creation of objects is profiled with the command
<code class="docutils notranslate"><span class="pre"><a href="../manualpages/Log/PetscLogObjectCreate.html">PetscLogObjectCreate</a>()</span></code></p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../manualpages/Log/PetscLogObjectCreate.html">PetscLogObjectCreate</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscObject.html">PetscObject</a></span><span class="w"> </span><span class="n">h</span><span class="p">);</span>
</pre></div>
</div>
<p>which logs the creation of any PETSc object. Just before an object is
destroyed, it should be logged with <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Log/PetscLogObjectDestroy.html">PetscLogObjectDestroy</a>()</span></code></p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../manualpages/Log/PetscLogObjectDestroy.html">PetscLogObjectDestroy</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscObject.html">PetscObject</a></span><span class="w"> </span><span class="n">h</span><span class="p">);</span>
</pre></div>
</div>
<p>These are called automatically by <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscHeaderCreate.html">PetscHeaderCreate</a>()</span></code> and
<code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscHeaderDestroy.html">PetscHeaderDestroy</a>()</span></code>, which are used in creating all objects
inherited from the basic object. Thus, these logging routines need never
be called directly.</p>
<p>It is also useful to log information about the state of an object, as
can be done with the command</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../manualpages/Log/PetscLogObjectState.html">PetscLogObjectState</a></span><span class="p">(</span><span class="n"><a href="../manualpages/Sys/PetscObject.html">PetscObject</a></span><span class="w"> </span><span class="n">h</span><span class="p">,</span><span class="k">const</span><span class="w"> </span><span class="kt">char</span><span class="w"> </span><span class="o">*</span><span class="n">format</span><span class="p">,...);</span>
</pre></div>
</div>
<p>For example, for sparse matrices we usually log the matrix dimensions
and number of nonzeros.</p>
</section>
<section id="profiling-events">
<h3>Profiling Events<a class="headerlink" href="#profiling-events" title="Link to this heading">#</a></h3>
<p>Events are logged by using the pair <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Log/PetscLogEventBegin.html">PetscLogEventBegin</a>()</span></code> and <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Log/PetscLogEventEnd.html">PetscLogEventEnd</a>()</span></code>.</p>
<p>This logging is usually done in the abstract interface file for the
operations, for example,
<a href="../src/mat/interface/matrix.c.html">src/mat/interface/matrix.c</a></p>
</section>
<section id="controlling-profiling">
<h3>Controlling Profiling<a class="headerlink" href="#controlling-profiling" title="Link to this heading">#</a></h3>
<p>Routines that control the default profiling available in PETSc include
the following</p>
<ul class="simple">
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Log/PetscLogDefaultBegin.html">PetscLogDefaultBegin</a>();</span></code></p></li>
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Log/PetscLogAllBegin.html">PetscLogAllBegin</a>();</span></code></p></li>
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Log/PetscLogDump.html">PetscLogDump</a>(const</span> <span class="pre">char</span> <span class="pre">*filename);</span></code></p></li>
<li><p><code class="docutils notranslate"><span class="pre"><a href="../manualpages/Log/PetscLogView.html">PetscLogView</a>(<a href="../manualpages/Viewer/PetscViewer.html">PetscViewer</a>);</span></code></p></li>
</ul>
<p>These routines are normally called by the <code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscInitialize.html">PetscInitialize</a>()</span></code> and
<code class="docutils notranslate"><span class="pre"><a href="../manualpages/Sys/PetscFinalize.html">PetscFinalize</a>()</span></code> routines when the option <code class="docutils notranslate"><span class="pre">-logview</span></code> is given.</p>
</section>
</section>
<section id="references">
<h2>References<a class="headerlink" href="#references" title="Link to this heading">#</a></h2>
<div class="docutils container" id="id4">
<div role="list" class="citation-list">
<div class="citation" id="id958" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id2">BBK+15</a><span class="fn-bracket">]</span></span>
<p>Satish Balay, Jed Brown, Matthew G. Knepley, Lois McInnes, and Barry Smith. Providing mixed language and legacy support within a library. In J. Carver, editor, <em>Software Engineering for Science</em>. Taylor &amp; Francis, 2015.</p>
</div>
<div class="citation" id="id821" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id1">BGMS98</a><span class="fn-bracket">]</span></span>
<p>Satish Balay, William D. Gropp, Lois Curfman McInnes, and Barry F. Smith. A microkernel design for component-based parallel numerical software systems. In <em>Proceedings of the SIAM Workshop on Object Oriented Methods for Inter-operable Scientific and Engineering Computing</em>, 58–67. SIAM, 1998. URL: <a class="reference external" href="ftp://info.mcs.anl.gov/pub/tech_reports/reports/P727.ps.Z">ftp://info.mcs.anl.gov/pub/tech_reports/reports/P727.ps.Z</a>.</p>
</div>
<div class="citation" id="id960" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id3">KBMS13</a><span class="fn-bracket">]</span></span>
<p>Matthew G. Knepley, Jed Brown, Lois Curfman McInnes, and Barry F. Smith. Accurately citing software and algorithms used in publications. In <em>First Workshop on Sustainable Software for Science: Practice and Experiences (WSSSPE), held at SC13</em>. 2013. <a class="reference external" href="https://doi.org/10.6084/m9.figshare.785731">doi:10.6084/m9.figshare.785731</a>.</p>
</div>
</div>
</div>
</section>
</section>


                </article>
              
              
              
              
              
                <footer class="prev-next-footer">
                  
<div class="prev-next-area">
    <a class="left-prev"
       href="design.html"
       title="previous page">
      <i class="fa-solid fa-angle-left"></i>
      <div class="prev-next-info">
        <p class="prev-next-subtitle">previous</p>
        <p class="prev-next-title">The Design of PETSc</p>
      </div>
    </a>
    <a class="right-next"
       href="objects.html"
       title="next page">
      <div class="prev-next-info">
        <p class="prev-next-subtitle">next</p>
        <p class="prev-next-title">Basic Object Design and Implementation</p>
      </div>
      <i class="fa-solid fa-angle-right"></i>
    </a>
</div>
                </footer>
              
            </div>
            
            
              
                <div class="bd-sidebar-secondary bd-toc"><div class="sidebar-secondary-items sidebar-secondary__inner">


  <div class="sidebar-secondary-item">
<div
    id="pst-page-navigation-heading-2"
    class="page-toc tocsection onthispage">
    <i class="fa-solid fa-list"></i> On this page
  </div>
  <nav class="bd-toc-nav page-toc" aria-labelledby="pst-page-navigation-heading-2">
    <ul class="visible nav section-nav flex-column">
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#petsc-types">PETSc Types</a></li>
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#implementation-of-error-handling">Implementation of Error Handling</a><ul class="nav section-nav flex-column">
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#simplified-interface">Simplified Interface</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#error-handlers">Error Handlers</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#error-codes">Error Codes</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#detailed-error-messages">Detailed Error Messages</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#c-exceptions">C++ Exceptions</a></li>
</ul>
</li>
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#memory-management">Memory Management</a></li>
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#implementation-of-profiling">Implementation of Profiling</a><ul class="nav section-nav flex-column">
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#profiling-object-creation-and-destruction">Profiling Object Creation and Destruction</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#profiling-events">Profiling Events</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#controlling-profiling">Controlling Profiling</a></li>
</ul>
</li>
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#references">References</a></li>
</ul>
  </nav></div>

  <div class="sidebar-secondary-item">

  
  <div class="tocsection editthispage">
    <a href="https://gitlab.com/petsc/petsc/-/edit/release/doc/developers/kernel.md">
      <i class="fa-solid fa-pencil"></i>
      
      
        
          Edit on GitLab
        
      
    </a>
  </div>
</div>

  <div class="sidebar-secondary-item">

  <div class="tocsection sourcelink">
    <a href="../_sources/developers/kernel.md.txt">
      <i class="fa-solid fa-file-lines"></i> Show Source
    </a>
  </div>
</div>

</div></div>
              
            
          </div>
          <footer class="bd-footer-content">
            
          </footer>
        
      </main>
    </div>
  </div>
  
  <!-- Scripts loaded after <body> so the DOM is not blocked -->
  <script src="../_static/scripts/bootstrap.js?digest=bd9e20870c6007c4c509"></script>
<script src="../_static/scripts/pydata-sphinx-theme.js?digest=bd9e20870c6007c4c509"></script>

  <footer class="bd-footer">
<div class="bd-footer__inner bd-page-width">
  
    <div class="footer-items__start">
      
        <div class="footer-item">

  <p class="copyright">
    
      © Copyright 1991-2025, UChicago Argonne, LLC and the PETSc Development Team.
      <br/>
    
  </p>
</div>
      
        <div class="footer-item">

  <p class="sphinx-version">
    Created using <a href="https://www.sphinx-doc.org/">Sphinx</a> 7.3.7.
    <br/>
  </p>
</div>
      
    </div>
  
  
  
    <div class="footer-items__end">
      
        <div class="footer-item">
<p class="theme-version">
  Built with the <a href="https://pydata-sphinx-theme.readthedocs.io/en/stable/index.html">PyData Sphinx Theme</a> 0.15.1.
</p></div>
      
        <div class="footer-item"><p class="last-updated">
  Last updated on 2025-04-30T13:10:40-0500 (v3.23.1).
  <br/>
</p></div>
      
    </div>
  
</div>

  </footer>
  </body>
</html>