File: petscmath.h.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 (1228 lines) | stat: -rw-r--r-- 125,824 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
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
<center><a href="https://gitlab.com/petsc/petsc/-/blob/966382dc56242773704ef5f5cee7aa2db3ebc577/include/petscmath.h">Actual source code: petscmath.h</a></center><br>

<html>
<head>
<title></title>
<meta name="generator" content="c2html 0.9.6">
<meta name="date" content="2025-04-30T18:14:50+00:00">
</head>

<body bgcolor="#FFFFFF">
<pre width=80>
<a name="line1">  1: </a><font color="#B22222">/*</font>
<a name="line2">  2: </a><font color="#B22222">    PETSc mathematics include file. Defines certain basic mathematical</font>
<a name="line3">  3: </a><font color="#B22222">    constants and functions for working with single, double, and quad precision</font>
<a name="line4">  4: </a><font color="#B22222">    floating point numbers as well as complex single and double.</font>

<a name="line6">  6: </a><font color="#B22222">    This file is included by petscsys.h and should not be used directly.</font>
<a name="line7">  7: </a><font color="#B22222">*/</font>
<a name="line8">  8: </a><font color="#A020F0">#pragma once</font>

<a name="line10"> 10: </a><font color="#A020F0">#include &lt;math.h&gt;</font>
<a name="line11"> 11: </a>#include <A href="../include/petscmacros.h.html">&lt;petscmacros.h&gt;</A>
<a name="line12"> 12: </a>#include <A href="../include/petscsystypes.h.html">&lt;petscsystypes.h&gt;</A>

<a name="line14"> 14: </a><font color="#B22222">/* SUBMANSEC = Sys */</font>

<a name="line16"> 16: </a><font color="#B22222">/*</font>
<a name="line17"> 17: </a><font color="#B22222">   Defines operations that are different for complex and real numbers.</font>
<a name="line18"> 18: </a><font color="#B22222">   All PETSc objects in one program are built around the object</font>
<a name="line19"> 19: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> which is either always a real or a complex.</font>
<a name="line20"> 20: </a><font color="#B22222">*/</font>

<a name="line22"> 22: </a><font color="#B22222">/*</font>
<a name="line23"> 23: </a><font color="#B22222">    Real number definitions</font>
<a name="line24"> 24: </a><font color="#B22222"> */</font>
<a name="line25"> 25: </a><font color="#A020F0">#if defined(PETSC_USE_REAL_SINGLE)</font>
<a name="line26"> 26: </a><strong><font color="#228B22">  #define PetscSqrtReal(a)        sqrtf(a)</font></strong>
<a name="line27"> 27: </a><strong><font color="#228B22">  #define PetscCbrtReal(a)        cbrtf(a)</font></strong>
<a name="line28"> 28: </a><strong><font color="#228B22">  #define PetscHypotReal(a, b)    hypotf(a, b)</font></strong>
<a name="line29"> 29: </a><strong><font color="#228B22">  #define PetscAtan2Real(a, b)    atan2f(a, b)</font></strong>
<a name="line30"> 30: </a><strong><font color="#228B22">  #define PetscPowReal(a, b)      powf(a, b)</font></strong>
<a name="line31"> 31: </a><strong><font color="#228B22">  #define PetscExpReal(a)         expf(a)</font></strong>
<a name="line32"> 32: </a><strong><font color="#228B22">  #define PetscLogReal(a)         logf(a)</font></strong>
<a name="line33"> 33: </a><strong><font color="#228B22">  #define PetscLog10Real(a)       log10f(a)</font></strong>
<a name="line34"> 34: </a><strong><font color="#228B22">  #define PetscLog2Real(a)        log2f(a)</font></strong>
<a name="line35"> 35: </a><strong><font color="#228B22">  #define PetscSinReal(a)         sinf(a)</font></strong>
<a name="line36"> 36: </a><strong><font color="#228B22">  #define PetscCosReal(a)         cosf(a)</font></strong>
<a name="line37"> 37: </a><strong><font color="#228B22">  #define PetscTanReal(a)         tanf(a)</font></strong>
<a name="line38"> 38: </a><strong><font color="#228B22">  #define PetscAsinReal(a)        asinf(a)</font></strong>
<a name="line39"> 39: </a><strong><font color="#228B22">  #define PetscAcosReal(a)        acosf(a)</font></strong>
<a name="line40"> 40: </a><strong><font color="#228B22">  #define PetscAtanReal(a)        atanf(a)</font></strong>
<a name="line41"> 41: </a><strong><font color="#228B22">  #define PetscSinhReal(a)        sinhf(a)</font></strong>
<a name="line42"> 42: </a><strong><font color="#228B22">  #define PetscCoshReal(a)        coshf(a)</font></strong>
<a name="line43"> 43: </a><strong><font color="#228B22">  #define PetscTanhReal(a)        tanhf(a)</font></strong>
<a name="line44"> 44: </a><strong><font color="#228B22">  #define PetscAsinhReal(a)       asinhf(a)</font></strong>
<a name="line45"> 45: </a><strong><font color="#228B22">  #define PetscAcoshReal(a)       acoshf(a)</font></strong>
<a name="line46"> 46: </a><strong><font color="#228B22">  #define PetscAtanhReal(a)       atanhf(a)</font></strong>
<a name="line47"> 47: </a><strong><font color="#228B22">  #define PetscErfReal(a)         erff(a)</font></strong>
<a name="line48"> 48: </a><strong><font color="#228B22">  #define PetscCeilReal(a)        ceilf(a)</font></strong>
<a name="line49"> 49: </a><strong><font color="#228B22">  #define PetscFloorReal(a)       floorf(a)</font></strong>
<a name="line50"> 50: </a><strong><font color="#228B22">  #define PetscRintReal(a)        rintf(a)</font></strong>
<a name="line51"> 51: </a><strong><font color="#228B22">  #define PetscFmodReal(a, b)     fmodf(a, b)</font></strong>
<a name="line52"> 52: </a><strong><font color="#228B22">  #define PetscCopysignReal(a, b) copysignf(a, b)</font></strong>
<a name="line53"> 53: </a><strong><font color="#228B22">  #define PetscTGamma(a)          tgammaf(a)</font></strong>
<a name="line54"> 54: </a><font color="#A020F0">  #if defined(PETSC_HAVE_LGAMMA_IS_GAMMA)</font>
<a name="line55"> 55: </a><strong><font color="#228B22">    #define PetscLGamma(a) gammaf(a)</font></strong>
<a name="line56"> 56: </a><font color="#A020F0">  #else</font>
<a name="line57"> 57: </a><strong><font color="#228B22">    #define PetscLGamma(a) lgammaf(a)</font></strong>
<a name="line58"> 58: </a><font color="#A020F0">  #endif</font>

<a name="line60"> 60: </a><font color="#A020F0">#elif defined(PETSC_USE_REAL_DOUBLE)</font>
<a name="line61"> 61: </a><strong><font color="#228B22">  #define PetscSqrtReal(a)        sqrt(a)</font></strong>
<a name="line62"> 62: </a><strong><font color="#228B22">  #define PetscCbrtReal(a)        cbrt(a)</font></strong>
<a name="line63"> 63: </a><strong><font color="#228B22">  #define PetscHypotReal(a, b)    hypot(a, b)</font></strong>
<a name="line64"> 64: </a><strong><font color="#228B22">  #define PetscAtan2Real(a, b)    atan2(a, b)</font></strong>
<a name="line65"> 65: </a><strong><font color="#228B22">  #define PetscPowReal(a, b)      pow(a, b)</font></strong>
<a name="line66"> 66: </a><strong><font color="#228B22">  #define PetscExpReal(a)         exp(a)</font></strong>
<a name="line67"> 67: </a><strong><font color="#228B22">  #define PetscLogReal(a)         log(a)</font></strong>
<a name="line68"> 68: </a><strong><font color="#228B22">  #define PetscLog10Real(a)       log10(a)</font></strong>
<a name="line69"> 69: </a><strong><font color="#228B22">  #define PetscLog2Real(a)        log2(a)</font></strong>
<a name="line70"> 70: </a><strong><font color="#228B22">  #define PetscSinReal(a)         sin(a)</font></strong>
<a name="line71"> 71: </a><strong><font color="#228B22">  #define PetscCosReal(a)         cos(a)</font></strong>
<a name="line72"> 72: </a><strong><font color="#228B22">  #define PetscTanReal(a)         tan(a)</font></strong>
<a name="line73"> 73: </a><strong><font color="#228B22">  #define PetscAsinReal(a)        asin(a)</font></strong>
<a name="line74"> 74: </a><strong><font color="#228B22">  #define PetscAcosReal(a)        acos(a)</font></strong>
<a name="line75"> 75: </a><strong><font color="#228B22">  #define PetscAtanReal(a)        atan(a)</font></strong>
<a name="line76"> 76: </a><strong><font color="#228B22">  #define PetscSinhReal(a)        sinh(a)</font></strong>
<a name="line77"> 77: </a><strong><font color="#228B22">  #define PetscCoshReal(a)        cosh(a)</font></strong>
<a name="line78"> 78: </a><strong><font color="#228B22">  #define PetscTanhReal(a)        tanh(a)</font></strong>
<a name="line79"> 79: </a><strong><font color="#228B22">  #define PetscAsinhReal(a)       asinh(a)</font></strong>
<a name="line80"> 80: </a><strong><font color="#228B22">  #define PetscAcoshReal(a)       acosh(a)</font></strong>
<a name="line81"> 81: </a><strong><font color="#228B22">  #define PetscAtanhReal(a)       atanh(a)</font></strong>
<a name="line82"> 82: </a><strong><font color="#228B22">  #define PetscErfReal(a)         erf(a)</font></strong>
<a name="line83"> 83: </a><strong><font color="#228B22">  #define PetscCeilReal(a)        ceil(a)</font></strong>
<a name="line84"> 84: </a><strong><font color="#228B22">  #define PetscFloorReal(a)       floor(a)</font></strong>
<a name="line85"> 85: </a><strong><font color="#228B22">  #define PetscRintReal(a)        rint(a)</font></strong>
<a name="line86"> 86: </a><strong><font color="#228B22">  #define PetscFmodReal(a, b)     fmod(a, b)</font></strong>
<a name="line87"> 87: </a><strong><font color="#228B22">  #define PetscCopysignReal(a, b) copysign(a, b)</font></strong>
<a name="line88"> 88: </a><strong><font color="#228B22">  #define PetscTGamma(a)          tgamma(a)</font></strong>
<a name="line89"> 89: </a><font color="#A020F0">  #if defined(PETSC_HAVE_LGAMMA_IS_GAMMA)</font>
<a name="line90"> 90: </a><strong><font color="#228B22">    #define PetscLGamma(a) gamma(a)</font></strong>
<a name="line91"> 91: </a><font color="#A020F0">  #else</font>
<a name="line92"> 92: </a><strong><font color="#228B22">    #define PetscLGamma(a) lgamma(a)</font></strong>
<a name="line93"> 93: </a><font color="#A020F0">  #endif</font>

<a name="line95"> 95: </a><font color="#A020F0">#elif defined(PETSC_USE_REAL___FLOAT128)</font>
<a name="line96"> 96: </a><strong><font color="#228B22">  #define PetscSqrtReal(a)        sqrtq(a)</font></strong>
<a name="line97"> 97: </a><strong><font color="#228B22">  #define PetscCbrtReal(a)        cbrtq(a)</font></strong>
<a name="line98"> 98: </a><strong><font color="#228B22">  #define PetscHypotReal(a, b)    hypotq(a, b)</font></strong>
<a name="line99"> 99: </a><strong><font color="#228B22">  #define PetscAtan2Real(a, b)    atan2q(a, b)</font></strong>
<a name="line100">100: </a><strong><font color="#228B22">  #define PetscPowReal(a, b)      powq(a, b)</font></strong>
<a name="line101">101: </a><strong><font color="#228B22">  #define PetscExpReal(a)         expq(a)</font></strong>
<a name="line102">102: </a><strong><font color="#228B22">  #define PetscLogReal(a)         logq(a)</font></strong>
<a name="line103">103: </a><strong><font color="#228B22">  #define PetscLog10Real(a)       log10q(a)</font></strong>
<a name="line104">104: </a><strong><font color="#228B22">  #define PetscLog2Real(a)        log2q(a)</font></strong>
<a name="line105">105: </a><strong><font color="#228B22">  #define PetscSinReal(a)         sinq(a)</font></strong>
<a name="line106">106: </a><strong><font color="#228B22">  #define PetscCosReal(a)         cosq(a)</font></strong>
<a name="line107">107: </a><strong><font color="#228B22">  #define PetscTanReal(a)         tanq(a)</font></strong>
<a name="line108">108: </a><strong><font color="#228B22">  #define PetscAsinReal(a)        asinq(a)</font></strong>
<a name="line109">109: </a><strong><font color="#228B22">  #define PetscAcosReal(a)        acosq(a)</font></strong>
<a name="line110">110: </a><strong><font color="#228B22">  #define PetscAtanReal(a)        atanq(a)</font></strong>
<a name="line111">111: </a><strong><font color="#228B22">  #define PetscSinhReal(a)        sinhq(a)</font></strong>
<a name="line112">112: </a><strong><font color="#228B22">  #define PetscCoshReal(a)        coshq(a)</font></strong>
<a name="line113">113: </a><strong><font color="#228B22">  #define PetscTanhReal(a)        tanhq(a)</font></strong>
<a name="line114">114: </a><strong><font color="#228B22">  #define PetscAsinhReal(a)       asinhq(a)</font></strong>
<a name="line115">115: </a><strong><font color="#228B22">  #define PetscAcoshReal(a)       acoshq(a)</font></strong>
<a name="line116">116: </a><strong><font color="#228B22">  #define PetscAtanhReal(a)       atanhq(a)</font></strong>
<a name="line117">117: </a><strong><font color="#228B22">  #define PetscErfReal(a)         erfq(a)</font></strong>
<a name="line118">118: </a><strong><font color="#228B22">  #define PetscCeilReal(a)        ceilq(a)</font></strong>
<a name="line119">119: </a><strong><font color="#228B22">  #define PetscFloorReal(a)       floorq(a)</font></strong>
<a name="line120">120: </a><strong><font color="#228B22">  #define PetscRintReal(a)        rintq(a)</font></strong>
<a name="line121">121: </a><strong><font color="#228B22">  #define PetscFmodReal(a, b)     fmodq(a, b)</font></strong>
<a name="line122">122: </a><strong><font color="#228B22">  #define PetscCopysignReal(a, b) copysignq(a, b)</font></strong>
<a name="line123">123: </a><strong><font color="#228B22">  #define PetscTGamma(a)          tgammaq(a)</font></strong>
<a name="line124">124: </a><font color="#A020F0">  #if defined(PETSC_HAVE_LGAMMA_IS_GAMMA)</font>
<a name="line125">125: </a><strong><font color="#228B22">    #define PetscLGamma(a) gammaq(a)</font></strong>
<a name="line126">126: </a><font color="#A020F0">  #else</font>
<a name="line127">127: </a><strong><font color="#228B22">    #define PetscLGamma(a) lgammaq(a)</font></strong>
<a name="line128">128: </a><font color="#A020F0">  #endif</font>

<a name="line130">130: </a><font color="#A020F0">#elif defined(PETSC_USE_REAL___FP16)</font>
<a name="line131">131: </a><strong><font color="#228B22">  #define PetscSqrtReal(a)        sqrtf(a)</font></strong>
<a name="line132">132: </a><strong><font color="#228B22">  #define PetscCbrtReal(a)        cbrtf(a)</font></strong>
<a name="line133">133: </a><strong><font color="#228B22">  #define PetscHypotReal(a, b)    hypotf(a, b)</font></strong>
<a name="line134">134: </a><strong><font color="#228B22">  #define PetscAtan2Real(a, b)    atan2f(a, b)</font></strong>
<a name="line135">135: </a><strong><font color="#228B22">  #define PetscPowReal(a, b)      powf(a, b)</font></strong>
<a name="line136">136: </a><strong><font color="#228B22">  #define PetscExpReal(a)         expf(a)</font></strong>
<a name="line137">137: </a><strong><font color="#228B22">  #define PetscLogReal(a)         logf(a)</font></strong>
<a name="line138">138: </a><strong><font color="#228B22">  #define PetscLog10Real(a)       log10f(a)</font></strong>
<a name="line139">139: </a><strong><font color="#228B22">  #define PetscLog2Real(a)        log2f(a)</font></strong>
<a name="line140">140: </a><strong><font color="#228B22">  #define PetscSinReal(a)         sinf(a)</font></strong>
<a name="line141">141: </a><strong><font color="#228B22">  #define PetscCosReal(a)         cosf(a)</font></strong>
<a name="line142">142: </a><strong><font color="#228B22">  #define PetscTanReal(a)         tanf(a)</font></strong>
<a name="line143">143: </a><strong><font color="#228B22">  #define PetscAsinReal(a)        asinf(a)</font></strong>
<a name="line144">144: </a><strong><font color="#228B22">  #define PetscAcosReal(a)        acosf(a)</font></strong>
<a name="line145">145: </a><strong><font color="#228B22">  #define PetscAtanReal(a)        atanf(a)</font></strong>
<a name="line146">146: </a><strong><font color="#228B22">  #define PetscSinhReal(a)        sinhf(a)</font></strong>
<a name="line147">147: </a><strong><font color="#228B22">  #define PetscCoshReal(a)        coshf(a)</font></strong>
<a name="line148">148: </a><strong><font color="#228B22">  #define PetscTanhReal(a)        tanhf(a)</font></strong>
<a name="line149">149: </a><strong><font color="#228B22">  #define PetscAsinhReal(a)       asinhf(a)</font></strong>
<a name="line150">150: </a><strong><font color="#228B22">  #define PetscAcoshReal(a)       acoshf(a)</font></strong>
<a name="line151">151: </a><strong><font color="#228B22">  #define PetscAtanhReal(a)       atanhf(a)</font></strong>
<a name="line152">152: </a><strong><font color="#228B22">  #define PetscErfReal(a)         erff(a)</font></strong>
<a name="line153">153: </a><strong><font color="#228B22">  #define PetscCeilReal(a)        ceilf(a)</font></strong>
<a name="line154">154: </a><strong><font color="#228B22">  #define PetscFloorReal(a)       floorf(a)</font></strong>
<a name="line155">155: </a><strong><font color="#228B22">  #define PetscRintReal(a)        rintf(a)</font></strong>
<a name="line156">156: </a><strong><font color="#228B22">  #define PetscFmodReal(a, b)     fmodf(a, b)</font></strong>
<a name="line157">157: </a><strong><font color="#228B22">  #define PetscCopysignReal(a, b) copysignf(a, b)</font></strong>
<a name="line158">158: </a><strong><font color="#228B22">  #define PetscTGamma(a)          tgammaf(a)</font></strong>
<a name="line159">159: </a><font color="#A020F0">  #if defined(PETSC_HAVE_LGAMMA_IS_GAMMA)</font>
<a name="line160">160: </a><strong><font color="#228B22">    #define PetscLGamma(a) gammaf(a)</font></strong>
<a name="line161">161: </a><font color="#A020F0">  #else</font>
<a name="line162">162: </a><strong><font color="#228B22">    #define PetscLGamma(a) lgammaf(a)</font></strong>
<a name="line163">163: </a><font color="#A020F0">  #endif</font>

<a name="line165">165: </a><font color="#A020F0">#endif </font><font color="#B22222">/* PETSC_USE_REAL_* */</font><font color="#A020F0"></font>

<a name="line167">167: </a><strong><font color="#4169E1"><a name="PetscSignReal"></a>static inline <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> PetscSignReal(<a href="../manualpages/Sys/PetscReal.html">PetscReal</a> a)</font></strong>
<a name="line168">168: </a>{
<a name="line169">169: </a>  <font color="#4169E1">return</font> (<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)((a &lt; (<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)0) ? -1 : ((a &gt; (<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)0) ? 1 : 0));
<a name="line170">170: </a>}

<a name="line172">172: </a><font color="#A020F0">#if !defined(PETSC_HAVE_LOG2)</font>
<a name="line173">173: </a><strong><font color="#228B22">  #undef PetscLog2Real</font></strong>
<a name="line174">174: </a><strong><font color="#4169E1"><a name="PetscLog2Real"></a>static inline <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> PetscLog2Real(<a href="../manualpages/Sys/PetscReal.html">PetscReal</a> a)</font></strong>
<a name="line175">175: </a>{
<a name="line176">176: </a>  <font color="#4169E1">return</font> PetscLogReal(a) / PetscLogReal((<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)2);
<a name="line177">177: </a>}
<a name="line178">178: </a><font color="#A020F0">#endif</font>

<a name="line180">180: </a><font color="#A020F0">#if defined(PETSC_HAVE_REAL___FLOAT128) &amp;&amp; !defined(PETSC_SKIP_REAL___FLOAT128)</font>
<a name="line181">181: </a><strong><font color="#4169E1">PETSC_EXTERN MPI_Datatype MPIU___FLOAT128 PETSC_ATTRIBUTE_MPI_TYPE_TAG(__float128)</font></strong>;
<a name="line182">182: </a><font color="#A020F0">#endif</font>
<a name="line183">183: </a><font color="#A020F0">#if defined(PETSC_HAVE_REAL___FP16) &amp;&amp; !defined(PETSC_SKIP_REAL___FP16)</font>
<a name="line184">184: </a><strong><font color="#4169E1">PETSC_EXTERN MPI_Datatype MPIU___FP16 PETSC_ATTRIBUTE_MPI_TYPE_TAG(__fp16)</font></strong>;
<a name="line185">185: </a><font color="#A020F0">#endif</font>

<a name="line187">187: </a><font color="#B22222">/*MC</font>
<a name="line188">188: </a><font color="#B22222">   <a href="../manualpages/Sys/MPIU_REAL.html">MPIU_REAL</a> - Portable MPI datatype corresponding to `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>` independent of what precision `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>` is in</font>

<a name="line190">190: </a><font color="#B22222">   Level: beginner</font>

<a name="line192">192: </a><font color="#B22222">   Note:</font>
<a name="line193">193: </a><font color="#B22222">   In MPI calls that require an MPI datatype that matches a `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>` or array of `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>` values, pass this value.</font>

<a name="line195">195: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>`, `<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a>`, `<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>`, `<a href="../manualpages/Sys/PetscInt.html">PetscInt</a>`, `<a href="../manualpages/Sys/MPIU_SCALAR.html">MPIU_SCALAR</a>`, `<a href="../manualpages/Sys/MPIU_COMPLEX.html">MPIU_COMPLEX</a>`, `<a href="../manualpages/Sys/MPIU_INT.html">MPIU_INT</a>`</font>
<a name="line196">196: </a><font color="#B22222">M*/</font>
<a name="line197">197: </a><font color="#A020F0">#if defined(PETSC_USE_REAL_SINGLE)</font>
<a name="line198">198: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/MPIU_REAL.html">MPIU_REAL</a> MPI_FLOAT</font></strong>
<a name="line199">199: </a><font color="#A020F0">#elif defined(PETSC_USE_REAL_DOUBLE)</font>
<a name="line200">200: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/MPIU_REAL.html">MPIU_REAL</a> MPI_DOUBLE</font></strong>
<a name="line201">201: </a><font color="#A020F0">#elif defined(PETSC_USE_REAL___FLOAT128)</font>
<a name="line202">202: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/MPIU_REAL.html">MPIU_REAL</a> MPIU___FLOAT128</font></strong>
<a name="line203">203: </a><font color="#A020F0">#elif defined(PETSC_USE_REAL___FP16)</font>
<a name="line204">204: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/MPIU_REAL.html">MPIU_REAL</a> MPIU___FP16</font></strong>
<a name="line205">205: </a><font color="#A020F0">#endif </font><font color="#B22222">/* PETSC_USE_REAL_* */</font><font color="#A020F0"></font>

<a name="line207">207: </a><font color="#B22222">/*</font>
<a name="line208">208: </a><font color="#B22222">    Complex number definitions</font>
<a name="line209">209: </a><font color="#B22222"> */</font>
<a name="line210">210: </a><font color="#A020F0">#if defined(PETSC_HAVE_COMPLEX)</font>
<a name="line211">211: </a><font color="#A020F0">  #if defined(__cplusplus) &amp;&amp; !defined(PETSC_USE_REAL___FLOAT128)</font>
<a name="line212">212: </a>  <font color="#B22222">/* C++ support of complex number */</font>

<a name="line214">214: </a><strong><font color="#228B22">    #define PetscRealPartComplex(a)      (static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a)).real()</font></strong>
<a name="line215">215: </a><strong><font color="#228B22">    #define PetscImaginaryPartComplex(a) (static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a)).imag()</font></strong>
<a name="line216">216: </a><strong><font color="#228B22">    #define PetscAbsComplex(a)           petsccomplexlib::abs(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>
<a name="line217">217: </a><strong><font color="#228B22">    #define PetscArgComplex(a)           petsccomplexlib::arg(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>
<a name="line218">218: </a><strong><font color="#228B22">    #define PetscConjComplex(a)          petsccomplexlib::conj(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>
<a name="line219">219: </a><strong><font color="#228B22">    #define PetscSqrtComplex(a)          petsccomplexlib::sqrt(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>
<a name="line220">220: </a><strong><font color="#228B22">    #define PetscPowComplex(a, b)        petsccomplexlib::pow(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a), static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(b))</font></strong>
<a name="line221">221: </a><strong><font color="#228B22">    #define PetscExpComplex(a)           petsccomplexlib::exp(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>
<a name="line222">222: </a><strong><font color="#228B22">    #define PetscLogComplex(a)           petsccomplexlib::log(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>
<a name="line223">223: </a><strong><font color="#228B22">    #define PetscSinComplex(a)           petsccomplexlib::sin(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>
<a name="line224">224: </a><strong><font color="#228B22">    #define PetscCosComplex(a)           petsccomplexlib::cos(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>
<a name="line225">225: </a><strong><font color="#228B22">    #define PetscTanComplex(a)           petsccomplexlib::tan(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>
<a name="line226">226: </a><strong><font color="#228B22">    #define PetscAsinComplex(a)          petsccomplexlib::asin(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>
<a name="line227">227: </a><strong><font color="#228B22">    #define PetscAcosComplex(a)          petsccomplexlib::acos(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>
<a name="line228">228: </a><strong><font color="#228B22">    #define PetscAtanComplex(a)          petsccomplexlib::atan(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>
<a name="line229">229: </a><strong><font color="#228B22">    #define PetscSinhComplex(a)          petsccomplexlib::sinh(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>
<a name="line230">230: </a><strong><font color="#228B22">    #define PetscCoshComplex(a)          petsccomplexlib::cosh(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>
<a name="line231">231: </a><strong><font color="#228B22">    #define PetscTanhComplex(a)          petsccomplexlib::tanh(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>
<a name="line232">232: </a><strong><font color="#228B22">    #define PetscAsinhComplex(a)         petsccomplexlib::asinh(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>
<a name="line233">233: </a><strong><font color="#228B22">    #define PetscAcoshComplex(a)         petsccomplexlib::acosh(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>
<a name="line234">234: </a><strong><font color="#228B22">    #define PetscAtanhComplex(a)         petsccomplexlib::atanh(static_cast&lt;<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>&gt;(a))</font></strong>

<a name="line236">236: </a>  <font color="#B22222">/* TODO: Add configure tests</font>

<a name="line238">238: </a><font color="#B22222">#if !defined(PETSC_HAVE_CXX_TAN_COMPLEX)</font>
<a name="line239">239: </a><font color="#B22222">#undef PetscTanComplex</font>
<a name="line240">240: </a><font color="#B22222">static inline <a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> PetscTanComplex(<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> z)</font>
<a name="line241">241: </a><font color="#B22222">{</font>
<a name="line242">242: </a><font color="#B22222">  return PetscSinComplex(z)/PetscCosComplex(z);</font>
<a name="line243">243: </a><font color="#B22222">}</font>
<a name="line244">244: </a><font color="#B22222">#endif</font>

<a name="line246">246: </a><font color="#B22222">#if !defined(PETSC_HAVE_CXX_TANH_COMPLEX)</font>
<a name="line247">247: </a><font color="#B22222">#undef PetscTanhComplex</font>
<a name="line248">248: </a><font color="#B22222">static inline <a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> PetscTanhComplex(<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> z)</font>
<a name="line249">249: </a><font color="#B22222">{</font>
<a name="line250">250: </a><font color="#B22222">  return PetscSinhComplex(z)/PetscCoshComplex(z);</font>
<a name="line251">251: </a><font color="#B22222">}</font>
<a name="line252">252: </a><font color="#B22222">#endif</font>

<a name="line254">254: </a><font color="#B22222">#if !defined(PETSC_HAVE_CXX_ASIN_COMPLEX)</font>
<a name="line255">255: </a><font color="#B22222">#undef PetscAsinComplex</font>
<a name="line256">256: </a><font color="#B22222">static inline <a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> PetscAsinComplex(<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> z)</font>
<a name="line257">257: </a><font color="#B22222">{</font>
<a name="line258">258: </a><font color="#B22222">  const <a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> j(0,1);</font>
<a name="line259">259: </a><font color="#B22222">  return -j*PetscLogComplex(j*z+PetscSqrtComplex(1.0f-z*z));</font>
<a name="line260">260: </a><font color="#B22222">}</font>
<a name="line261">261: </a><font color="#B22222">#endif</font>

<a name="line263">263: </a><font color="#B22222">#if !defined(PETSC_HAVE_CXX_ACOS_COMPLEX)</font>
<a name="line264">264: </a><font color="#B22222">#undef PetscAcosComplex</font>
<a name="line265">265: </a><font color="#B22222">static inline <a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> PetscAcosComplex(<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> z)</font>
<a name="line266">266: </a><font color="#B22222">{</font>
<a name="line267">267: </a><font color="#B22222">  const <a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> j(0,1);</font>
<a name="line268">268: </a><font color="#B22222">  return j*PetscLogComplex(z-j*PetscSqrtComplex(1.0f-z*z));</font>
<a name="line269">269: </a><font color="#B22222">}</font>
<a name="line270">270: </a><font color="#B22222">#endif</font>

<a name="line272">272: </a><font color="#B22222">#if !defined(PETSC_HAVE_CXX_ATAN_COMPLEX)</font>
<a name="line273">273: </a><font color="#B22222">#undef PetscAtanComplex</font>
<a name="line274">274: </a><font color="#B22222">static inline <a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> PetscAtanComplex(<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> z)</font>
<a name="line275">275: </a><font color="#B22222">{</font>
<a name="line276">276: </a><font color="#B22222">  const <a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> j(0,1);</font>
<a name="line277">277: </a><font color="#B22222">  return 0.5f*j*PetscLogComplex((1.0f-j*z)/(1.0f+j*z));</font>
<a name="line278">278: </a><font color="#B22222">}</font>
<a name="line279">279: </a><font color="#B22222">#endif</font>

<a name="line281">281: </a><font color="#B22222">#if !defined(PETSC_HAVE_CXX_ASINH_COMPLEX)</font>
<a name="line282">282: </a><font color="#B22222">#undef PetscAsinhComplex</font>
<a name="line283">283: </a><font color="#B22222">static inline <a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> PetscAsinhComplex(<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> z)</font>
<a name="line284">284: </a><font color="#B22222">{</font>
<a name="line285">285: </a><font color="#B22222">  return PetscLogComplex(z+PetscSqrtComplex(z*z+1.0f));</font>
<a name="line286">286: </a><font color="#B22222">}</font>
<a name="line287">287: </a><font color="#B22222">#endif</font>

<a name="line289">289: </a><font color="#B22222">#if !defined(PETSC_HAVE_CXX_ACOSH_COMPLEX)</font>
<a name="line290">290: </a><font color="#B22222">#undef PetscAcoshComplex</font>
<a name="line291">291: </a><font color="#B22222">static inline <a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> PetscAcoshComplex(<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> z)</font>
<a name="line292">292: </a><font color="#B22222">{</font>
<a name="line293">293: </a><font color="#B22222">  return PetscLogComplex(z+PetscSqrtComplex(z*z-1.0f));</font>
<a name="line294">294: </a><font color="#B22222">}</font>
<a name="line295">295: </a><font color="#B22222">#endif</font>

<a name="line297">297: </a><font color="#B22222">#if !defined(PETSC_HAVE_CXX_ATANH_COMPLEX)</font>
<a name="line298">298: </a><font color="#B22222">#undef PetscAtanhComplex</font>
<a name="line299">299: </a><font color="#B22222">static inline <a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> PetscAtanhComplex(<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> z)</font>
<a name="line300">300: </a><font color="#B22222">{</font>
<a name="line301">301: </a><font color="#B22222">  return 0.5f*PetscLogComplex((1.0f+z)/(1.0f-z));</font>
<a name="line302">302: </a><font color="#B22222">}</font>
<a name="line303">303: </a><font color="#B22222">#endif</font>

<a name="line305">305: </a><font color="#B22222">*/</font>

<a name="line307">307: </a><font color="#A020F0">  #else </font><font color="#B22222">/* C99 support of complex number */</font><font color="#A020F0"></font>

<a name="line309">309: </a><font color="#A020F0">    #if defined(PETSC_USE_REAL_SINGLE)</font>
<a name="line310">310: </a><strong><font color="#228B22">      #define PetscRealPartComplex(a)      crealf(a)</font></strong>
<a name="line311">311: </a><strong><font color="#228B22">      #define PetscImaginaryPartComplex(a) cimagf(a)</font></strong>
<a name="line312">312: </a><strong><font color="#228B22">      #define PetscAbsComplex(a)           cabsf(a)</font></strong>
<a name="line313">313: </a><strong><font color="#228B22">      #define PetscArgComplex(a)           cargf(a)</font></strong>
<a name="line314">314: </a><strong><font color="#228B22">      #define PetscConjComplex(a)          conjf(a)</font></strong>
<a name="line315">315: </a><strong><font color="#228B22">      #define PetscSqrtComplex(a)          csqrtf(a)</font></strong>
<a name="line316">316: </a><strong><font color="#228B22">      #define PetscPowComplex(a, b)        cpowf(a, b)</font></strong>
<a name="line317">317: </a><strong><font color="#228B22">      #define PetscExpComplex(a)           cexpf(a)</font></strong>
<a name="line318">318: </a><strong><font color="#228B22">      #define PetscLogComplex(a)           clogf(a)</font></strong>
<a name="line319">319: </a><strong><font color="#228B22">      #define PetscSinComplex(a)           csinf(a)</font></strong>
<a name="line320">320: </a><strong><font color="#228B22">      #define PetscCosComplex(a)           ccosf(a)</font></strong>
<a name="line321">321: </a><strong><font color="#228B22">      #define PetscTanComplex(a)           ctanf(a)</font></strong>
<a name="line322">322: </a><strong><font color="#228B22">      #define PetscAsinComplex(a)          casinf(a)</font></strong>
<a name="line323">323: </a><strong><font color="#228B22">      #define PetscAcosComplex(a)          cacosf(a)</font></strong>
<a name="line324">324: </a><strong><font color="#228B22">      #define PetscAtanComplex(a)          catanf(a)</font></strong>
<a name="line325">325: </a><strong><font color="#228B22">      #define PetscSinhComplex(a)          csinhf(a)</font></strong>
<a name="line326">326: </a><strong><font color="#228B22">      #define PetscCoshComplex(a)          ccoshf(a)</font></strong>
<a name="line327">327: </a><strong><font color="#228B22">      #define PetscTanhComplex(a)          ctanhf(a)</font></strong>
<a name="line328">328: </a><strong><font color="#228B22">      #define PetscAsinhComplex(a)         casinhf(a)</font></strong>
<a name="line329">329: </a><strong><font color="#228B22">      #define PetscAcoshComplex(a)         cacoshf(a)</font></strong>
<a name="line330">330: </a><strong><font color="#228B22">      #define PetscAtanhComplex(a)         catanhf(a)</font></strong>

<a name="line332">332: </a><font color="#A020F0">    #elif defined(PETSC_USE_REAL_DOUBLE)</font>
<a name="line333">333: </a><strong><font color="#228B22">      #define PetscRealPartComplex(a)      creal(a)</font></strong>
<a name="line334">334: </a><strong><font color="#228B22">      #define PetscImaginaryPartComplex(a) cimag(a)</font></strong>
<a name="line335">335: </a><strong><font color="#228B22">      #define PetscAbsComplex(a)           cabs(a)</font></strong>
<a name="line336">336: </a><strong><font color="#228B22">      #define PetscArgComplex(a)           carg(a)</font></strong>
<a name="line337">337: </a><strong><font color="#228B22">      #define PetscConjComplex(a)          conj(a)</font></strong>
<a name="line338">338: </a><strong><font color="#228B22">      #define PetscSqrtComplex(a)          csqrt(a)</font></strong>
<a name="line339">339: </a><strong><font color="#228B22">      #define PetscPowComplex(a, b)        cpow(a, b)</font></strong>
<a name="line340">340: </a><strong><font color="#228B22">      #define PetscExpComplex(a)           cexp(a)</font></strong>
<a name="line341">341: </a><strong><font color="#228B22">      #define PetscLogComplex(a)           clog(a)</font></strong>
<a name="line342">342: </a><strong><font color="#228B22">      #define PetscSinComplex(a)           csin(a)</font></strong>
<a name="line343">343: </a><strong><font color="#228B22">      #define PetscCosComplex(a)           ccos(a)</font></strong>
<a name="line344">344: </a><strong><font color="#228B22">      #define PetscTanComplex(a)           ctan(a)</font></strong>
<a name="line345">345: </a><strong><font color="#228B22">      #define PetscAsinComplex(a)          casin(a)</font></strong>
<a name="line346">346: </a><strong><font color="#228B22">      #define PetscAcosComplex(a)          cacos(a)</font></strong>
<a name="line347">347: </a><strong><font color="#228B22">      #define PetscAtanComplex(a)          catan(a)</font></strong>
<a name="line348">348: </a><strong><font color="#228B22">      #define PetscSinhComplex(a)          csinh(a)</font></strong>
<a name="line349">349: </a><strong><font color="#228B22">      #define PetscCoshComplex(a)          ccosh(a)</font></strong>
<a name="line350">350: </a><strong><font color="#228B22">      #define PetscTanhComplex(a)          ctanh(a)</font></strong>
<a name="line351">351: </a><strong><font color="#228B22">      #define PetscAsinhComplex(a)         casinh(a)</font></strong>
<a name="line352">352: </a><strong><font color="#228B22">      #define PetscAcoshComplex(a)         cacosh(a)</font></strong>
<a name="line353">353: </a><strong><font color="#228B22">      #define PetscAtanhComplex(a)         catanh(a)</font></strong>

<a name="line355">355: </a><font color="#A020F0">    #elif defined(PETSC_USE_REAL___FLOAT128)</font>
<a name="line356">356: </a><strong><font color="#228B22">      #define PetscRealPartComplex(a)      crealq(a)</font></strong>
<a name="line357">357: </a><strong><font color="#228B22">      #define PetscImaginaryPartComplex(a) cimagq(a)</font></strong>
<a name="line358">358: </a><strong><font color="#228B22">      #define PetscAbsComplex(a)           cabsq(a)</font></strong>
<a name="line359">359: </a><strong><font color="#228B22">      #define PetscArgComplex(a)           cargq(a)</font></strong>
<a name="line360">360: </a><strong><font color="#228B22">      #define PetscConjComplex(a)          conjq(a)</font></strong>
<a name="line361">361: </a><strong><font color="#228B22">      #define PetscSqrtComplex(a)          csqrtq(a)</font></strong>
<a name="line362">362: </a><strong><font color="#228B22">      #define PetscPowComplex(a, b)        cpowq(a, b)</font></strong>
<a name="line363">363: </a><strong><font color="#228B22">      #define PetscExpComplex(a)           cexpq(a)</font></strong>
<a name="line364">364: </a><strong><font color="#228B22">      #define PetscLogComplex(a)           clogq(a)</font></strong>
<a name="line365">365: </a><strong><font color="#228B22">      #define PetscSinComplex(a)           csinq(a)</font></strong>
<a name="line366">366: </a><strong><font color="#228B22">      #define PetscCosComplex(a)           ccosq(a)</font></strong>
<a name="line367">367: </a><strong><font color="#228B22">      #define PetscTanComplex(a)           ctanq(a)</font></strong>
<a name="line368">368: </a><strong><font color="#228B22">      #define PetscAsinComplex(a)          casinq(a)</font></strong>
<a name="line369">369: </a><strong><font color="#228B22">      #define PetscAcosComplex(a)          cacosq(a)</font></strong>
<a name="line370">370: </a><strong><font color="#228B22">      #define PetscAtanComplex(a)          catanq(a)</font></strong>
<a name="line371">371: </a><strong><font color="#228B22">      #define PetscSinhComplex(a)          csinhq(a)</font></strong>
<a name="line372">372: </a><strong><font color="#228B22">      #define PetscCoshComplex(a)          ccoshq(a)</font></strong>
<a name="line373">373: </a><strong><font color="#228B22">      #define PetscTanhComplex(a)          ctanhq(a)</font></strong>
<a name="line374">374: </a><strong><font color="#228B22">      #define PetscAsinhComplex(a)         casinhq(a)</font></strong>
<a name="line375">375: </a><strong><font color="#228B22">      #define PetscAcoshComplex(a)         cacoshq(a)</font></strong>
<a name="line376">376: </a><strong><font color="#228B22">      #define PetscAtanhComplex(a)         catanhq(a)</font></strong>

<a name="line378">378: </a><font color="#A020F0">    #endif </font><font color="#B22222">/* PETSC_USE_REAL_* */</font><font color="#A020F0"></font>
<a name="line379">379: </a><font color="#A020F0">  #endif   </font><font color="#B22222">/* (__cplusplus) */</font><font color="#A020F0"></font>

<a name="line381">381: </a><font color="#B22222">/*MC</font>
<a name="line382">382: </a><font color="#B22222">   <a href="../manualpages/Sys/PETSC_i.html">PETSC_i</a> - the pure imaginary complex number i</font>

<a name="line384">384: </a><font color="#B22222">   Level: intermediate</font>

<a name="line386">386: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>`, `<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a>`</font>
<a name="line387">387: </a><font color="#B22222">M*/</font>
<a name="line388">388: </a>PETSC_EXTERN <a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> <a href="../manualpages/Sys/PETSC_i.html">PETSC_i</a>;

<a name="line390">390: </a><font color="#B22222">/*</font>
<a name="line391">391: </a><font color="#B22222">   Try to do the right thing for complex number construction: see</font>
<a name="line392">392: </a><font color="#B22222">   http://www.open-std.org/jtc1/sc22/wg14/www/docs/n1464.htm</font>
<a name="line393">393: </a><font color="#B22222">   for details</font>
<a name="line394">394: </a><font color="#B22222">*/</font>
<a name="line395">395: </a><strong><font color="#4169E1"><a name="PetscCMPLX"></a>static inline <a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> PetscCMPLX(<a href="../manualpages/Sys/PetscReal.html">PetscReal</a> x, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> y)</font></strong>
<a name="line396">396: </a>{
<a name="line397">397: </a><font color="#A020F0">  #if defined(__cplusplus) &amp;&amp; !defined(PETSC_USE_REAL___FLOAT128)</font>
<a name="line398">398: </a>  <font color="#4169E1">return</font> <a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>(x, y);
<a name="line399">399: </a><font color="#A020F0">  #elif defined(_Imaginary_I)</font>
<a name="line400">400: </a>  <font color="#4169E1">return</font> x + y * _Imaginary_I;
<a name="line401">401: </a><font color="#A020F0">  #else</font>
<a name="line402">402: </a>  { <font color="#B22222">/* In both C99 and C11 (ISO/IEC 9899, Section 6.2.5),</font>

<a name="line404">404: </a><font color="#B22222">       "For each floating type there is a corresponding real type, which is always a real floating</font>
<a name="line405">405: </a><font color="#B22222">       type. For real floating types, it is the same type. For complex types, it is the type given</font>
<a name="line406">406: </a><font color="#B22222">       by deleting the keyword _Complex from the type name."</font>

<a name="line408">408: </a><font color="#B22222">       So type punning should be portable. */</font>
<a name="line409">409: </a>    <font color="#4169E1">union</font>
<a name="line410">410: </a>    {
<a name="line411">411: </a>      <a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a> z;
<a name="line412">412: </a>      <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>    f[2];
<a name="line413">413: </a>    } uz;

<a name="line415">415: </a>    uz.f[0] = x;
<a name="line416">416: </a>    uz.f[1] = y;
<a name="line417">417: </a>    <font color="#4169E1">return</font> uz.z;
<a name="line418">418: </a>  }
<a name="line419">419: </a><font color="#A020F0">  #endif</font>
<a name="line420">420: </a>}

<a name="line422">422: </a><strong><font color="#228B22">  #define MPIU_C_COMPLEX        MPI_C_COMPLEX PETSC_DEPRECATED_MACRO(3, 15, 0, </font><font color="#666666">"MPI_C_COMPLEX"</font><font color="#228B22">, )</font></strong>
<a name="line423">423: </a><strong><font color="#228B22">  #define MPIU_C_DOUBLE_COMPLEX MPI_C_DOUBLE_COMPLEX PETSC_DEPRECATED_MACRO(3, 15, 0, </font><font color="#666666">"MPI_C_DOUBLE_COMPLEX"</font><font color="#228B22">, )</font></strong>

<a name="line425">425: </a><font color="#A020F0">  #if defined(PETSC_HAVE_REAL___FLOAT128) &amp;&amp; !defined(PETSC_SKIP_REAL___FLOAT128)</font>
<a name="line426">426: </a>    <font color="#B22222">// if complex is not used, then quadmath.h won't be included by petscsystypes.h</font>
<a name="line427">427: </a><font color="#A020F0">    #if defined(PETSC_USE_COMPLEX)</font>
<a name="line428">428: </a><strong><font color="#228B22">      #define MPIU___COMPLEX128_ATTR_TAG PETSC_ATTRIBUTE_MPI_TYPE_TAG(__complex128)</font></strong>
<a name="line429">429: </a><font color="#A020F0">    #else</font>
<a name="line430">430: </a><strong><font color="#228B22">      #define MPIU___COMPLEX128_ATTR_TAG</font></strong>
<a name="line431">431: </a><font color="#A020F0">    #endif</font>

<a name="line433">433: </a>PETSC_EXTERN MPI_Datatype MPIU___COMPLEX128 MPIU___COMPLEX128_ATTR_TAG;

<a name="line435">435: </a><strong><font color="#228B22">    #undef MPIU___COMPLEX128_ATTR_TAG</font></strong>
<a name="line436">436: </a><font color="#A020F0">  #endif </font><font color="#B22222">/* PETSC_HAVE_REAL___FLOAT128 */</font><font color="#A020F0"></font>

<a name="line438">438: </a>  <font color="#B22222">/*MC</font>
<a name="line439">439: </a><font color="#B22222">   <a href="../manualpages/Sys/MPIU_COMPLEX.html">MPIU_COMPLEX</a> - Portable MPI datatype corresponding to `<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>` independent of the precision of `<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>`</font>

<a name="line441">441: </a><font color="#B22222">   Level: beginner</font>

<a name="line443">443: </a><font color="#B22222">   Note:</font>
<a name="line444">444: </a><font color="#B22222">   In MPI calls that require an MPI datatype that matches a `<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>` or array of `<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>` values, pass this value.</font>

<a name="line446">446: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>`, `<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a>`, `<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>`, `<a href="../manualpages/Sys/PetscInt.html">PetscInt</a>`, `<a href="../manualpages/Sys/MPIU_REAL.html">MPIU_REAL</a>`, `<a href="../manualpages/Sys/MPIU_SCALAR.html">MPIU_SCALAR</a>`, `<a href="../manualpages/Sys/MPIU_COMPLEX.html">MPIU_COMPLEX</a>`, `<a href="../manualpages/Sys/MPIU_INT.html">MPIU_INT</a>`, `<a href="../manualpages/Sys/PETSC_i.html">PETSC_i</a>`</font>
<a name="line447">447: </a><font color="#B22222">M*/</font>
<a name="line448">448: </a><font color="#A020F0">  #if defined(PETSC_USE_REAL_SINGLE)</font>
<a name="line449">449: </a><strong><font color="#228B22">    #define <a href="../manualpages/Sys/MPIU_COMPLEX.html">MPIU_COMPLEX</a> MPI_C_COMPLEX</font></strong>
<a name="line450">450: </a><font color="#A020F0">  #elif defined(PETSC_USE_REAL_DOUBLE)</font>
<a name="line451">451: </a><strong><font color="#228B22">    #define <a href="../manualpages/Sys/MPIU_COMPLEX.html">MPIU_COMPLEX</a> MPI_C_DOUBLE_COMPLEX</font></strong>
<a name="line452">452: </a><font color="#A020F0">  #elif defined(PETSC_USE_REAL___FLOAT128)</font>
<a name="line453">453: </a><strong><font color="#228B22">    #define <a href="../manualpages/Sys/MPIU_COMPLEX.html">MPIU_COMPLEX</a> MPIU___COMPLEX128</font></strong>
<a name="line454">454: </a><font color="#A020F0">  #elif defined(PETSC_USE_REAL___FP16)</font>
<a name="line455">455: </a><strong><font color="#228B22">    #define <a href="../manualpages/Sys/MPIU_COMPLEX.html">MPIU_COMPLEX</a> MPI_C_COMPLEX</font></strong>
<a name="line456">456: </a><font color="#A020F0">  #endif </font><font color="#B22222">/* PETSC_USE_REAL_* */</font><font color="#A020F0"></font>

<a name="line458">458: </a><font color="#A020F0">#endif </font><font color="#B22222">/* PETSC_HAVE_COMPLEX */</font><font color="#A020F0"></font>

<a name="line460">460: </a><font color="#B22222">/*</font>
<a name="line461">461: </a><font color="#B22222">    Scalar number definitions</font>
<a name="line462">462: </a><font color="#B22222"> */</font>
<a name="line463">463: </a><font color="#A020F0">#if defined(PETSC_USE_COMPLEX) &amp;&amp; defined(PETSC_HAVE_COMPLEX)</font>
<a name="line464">464: </a>  <font color="#B22222">/*MC</font>
<a name="line465">465: </a><font color="#B22222">   <a href="../manualpages/Sys/MPIU_SCALAR.html">MPIU_SCALAR</a> - Portable MPI datatype corresponding to `<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a>` independent of the precision of `<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a>`</font>

<a name="line467">467: </a><font color="#B22222">   Level: beginner</font>

<a name="line469">469: </a><font color="#B22222">   Note:</font>
<a name="line470">470: </a><font color="#B22222">   In MPI calls that require an MPI datatype that matches a `<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a>` or array of `<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a>` values, pass this value.</font>

<a name="line472">472: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>`, `<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a>`, `<a href="../manualpages/Sys/PetscComplex.html">PetscComplex</a>`, `<a href="../manualpages/Sys/PetscInt.html">PetscInt</a>`, `<a href="../manualpages/Sys/MPIU_REAL.html">MPIU_REAL</a>`, `<a href="../manualpages/Sys/MPIU_COMPLEX.html">MPIU_COMPLEX</a>`, `<a href="../manualpages/Sys/MPIU_INT.html">MPIU_INT</a>`</font>
<a name="line473">473: </a><font color="#B22222">M*/</font>
<a name="line474">474: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/MPIU_SCALAR.html">MPIU_SCALAR</a> <a href="../manualpages/Sys/MPIU_COMPLEX.html">MPIU_COMPLEX</a></font></strong>

<a name="line476">476: </a>  <font color="#B22222">/*MC</font>
<a name="line477">477: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscRealPart.html">PetscRealPart</a> - Returns the real part of a `<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a>`</font>

<a name="line479">479: </a><font color="#B22222">   Synopsis:</font>
<a name="line480">480: </a>#include <A href="../include/petscmath.h.html">&lt;petscmath.h&gt;</A>
<a name="line481">481: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> <a href="../manualpages/Sys/PetscRealPart.html">PetscRealPart</a>(<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> v)</font>

<a name="line483">483: </a><font color="#B22222">   Not Collective</font>

<a name="line485">485: </a><font color="#B22222">   Input Parameter:</font>
<a name="line486">486: </a><font color="#B22222">.  v - value to find the real part of</font>

<a name="line488">488: </a><font color="#B22222">   Level: beginner</font>

<a name="line490">490: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a>`, `<a href="../manualpages/Sys/PetscImaginaryPart.html">PetscImaginaryPart</a>()`, `<a href="../manualpages/Sys/PetscMax.html">PetscMax</a>()`, `<a href="../manualpages/Sys/PetscClipInterval.html">PetscClipInterval</a>()`, `<a href="../manualpages/Sys/PetscAbsInt.html">PetscAbsInt</a>()`, `<a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>()`, `<a href="../manualpages/Sys/PetscSqr.html">PetscSqr</a>()`</font>
<a name="line491">491: </a><font color="#B22222">M*/</font>
<a name="line492">492: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PetscRealPart.html">PetscRealPart</a>(a) PetscRealPartComplex(a)</font></strong>

<a name="line494">494: </a>  <font color="#B22222">/*MC</font>
<a name="line495">495: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscImaginaryPart.html">PetscImaginaryPart</a> - Returns the imaginary part of a `<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a>`</font>

<a name="line497">497: </a><font color="#B22222">   Synopsis:</font>
<a name="line498">498: </a>#include <A href="../include/petscmath.h.html">&lt;petscmath.h&gt;</A>
<a name="line499">499: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> <a href="../manualpages/Sys/PetscImaginaryPart.html">PetscImaginaryPart</a>(<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> v)</font>

<a name="line501">501: </a><font color="#B22222">   Not Collective</font>

<a name="line503">503: </a><font color="#B22222">   Input Parameter:</font>
<a name="line504">504: </a><font color="#B22222">.  v - value to find the imaginary part of</font>

<a name="line506">506: </a><font color="#B22222">   Level: beginner</font>

<a name="line508">508: </a><font color="#B22222">   Note:</font>
<a name="line509">509: </a><font color="#B22222">   If PETSc was configured for real numbers then this always returns the value 0</font>

<a name="line511">511: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a>`, `<a href="../manualpages/Sys/PetscRealPart.html">PetscRealPart</a>()`, `<a href="../manualpages/Sys/PetscMax.html">PetscMax</a>()`, `<a href="../manualpages/Sys/PetscClipInterval.html">PetscClipInterval</a>()`, `<a href="../manualpages/Sys/PetscAbsInt.html">PetscAbsInt</a>()`, `<a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>()`, `<a href="../manualpages/Sys/PetscSqr.html">PetscSqr</a>()`</font>
<a name="line512">512: </a><font color="#B22222">M*/</font>
<a name="line513">513: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PetscImaginaryPart.html">PetscImaginaryPart</a>(a) PetscImaginaryPartComplex(a)</font></strong>

<a name="line515">515: </a><strong><font color="#228B22">  #define PetscAbsScalar(a)    PetscAbsComplex(a)</font></strong>
<a name="line516">516: </a><strong><font color="#228B22">  #define PetscArgScalar(a)    PetscArgComplex(a)</font></strong>
<a name="line517">517: </a><strong><font color="#228B22">  #define PetscConj(a)         PetscConjComplex(a)</font></strong>
<a name="line518">518: </a><strong><font color="#228B22">  #define PetscSqrtScalar(a)   PetscSqrtComplex(a)</font></strong>
<a name="line519">519: </a><strong><font color="#228B22">  #define PetscPowScalar(a, b) PetscPowComplex(a, b)</font></strong>
<a name="line520">520: </a><strong><font color="#228B22">  #define PetscExpScalar(a)    PetscExpComplex(a)</font></strong>
<a name="line521">521: </a><strong><font color="#228B22">  #define PetscLogScalar(a)    PetscLogComplex(a)</font></strong>
<a name="line522">522: </a><strong><font color="#228B22">  #define PetscSinScalar(a)    PetscSinComplex(a)</font></strong>
<a name="line523">523: </a><strong><font color="#228B22">  #define PetscCosScalar(a)    PetscCosComplex(a)</font></strong>
<a name="line524">524: </a><strong><font color="#228B22">  #define PetscTanScalar(a)    PetscTanComplex(a)</font></strong>
<a name="line525">525: </a><strong><font color="#228B22">  #define PetscAsinScalar(a)   PetscAsinComplex(a)</font></strong>
<a name="line526">526: </a><strong><font color="#228B22">  #define PetscAcosScalar(a)   PetscAcosComplex(a)</font></strong>
<a name="line527">527: </a><strong><font color="#228B22">  #define PetscAtanScalar(a)   PetscAtanComplex(a)</font></strong>
<a name="line528">528: </a><strong><font color="#228B22">  #define PetscSinhScalar(a)   PetscSinhComplex(a)</font></strong>
<a name="line529">529: </a><strong><font color="#228B22">  #define PetscCoshScalar(a)   PetscCoshComplex(a)</font></strong>
<a name="line530">530: </a><strong><font color="#228B22">  #define PetscTanhScalar(a)   PetscTanhComplex(a)</font></strong>
<a name="line531">531: </a><strong><font color="#228B22">  #define PetscAsinhScalar(a)  PetscAsinhComplex(a)</font></strong>
<a name="line532">532: </a><strong><font color="#228B22">  #define PetscAcoshScalar(a)  PetscAcoshComplex(a)</font></strong>
<a name="line533">533: </a><strong><font color="#228B22">  #define PetscAtanhScalar(a)  PetscAtanhComplex(a)</font></strong>

<a name="line535">535: </a><font color="#A020F0">#else </font><font color="#B22222">/* PETSC_USE_COMPLEX */</font><font color="#A020F0"></font>
<a name="line536">536: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/MPIU_SCALAR.html">MPIU_SCALAR</a>           <a href="../manualpages/Sys/MPIU_REAL.html">MPIU_REAL</a></font></strong>
<a name="line537">537: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PetscRealPart.html">PetscRealPart</a>(a)      (a)</font></strong>
<a name="line538">538: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PetscImaginaryPart.html">PetscImaginaryPart</a>(a) ((<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)0)</font></strong>
<a name="line539">539: </a><strong><font color="#228B22">  #define PetscAbsScalar(a)     <a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>(a)</font></strong>
<a name="line540">540: </a><strong><font color="#228B22">  #define PetscArgScalar(a)     (((a) &lt; (<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)0) ? <a href="../manualpages/Sys/PETSC_PI.html">PETSC_PI</a> : (<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)0)</font></strong>
<a name="line541">541: </a><strong><font color="#228B22">  #define PetscConj(a)          (a)</font></strong>
<a name="line542">542: </a><strong><font color="#228B22">  #define PetscSqrtScalar(a)    PetscSqrtReal(a)</font></strong>
<a name="line543">543: </a><strong><font color="#228B22">  #define PetscPowScalar(a, b)  PetscPowReal(a, b)</font></strong>
<a name="line544">544: </a><strong><font color="#228B22">  #define PetscExpScalar(a)     PetscExpReal(a)</font></strong>
<a name="line545">545: </a><strong><font color="#228B22">  #define PetscLogScalar(a)     PetscLogReal(a)</font></strong>
<a name="line546">546: </a><strong><font color="#228B22">  #define PetscSinScalar(a)     PetscSinReal(a)</font></strong>
<a name="line547">547: </a><strong><font color="#228B22">  #define PetscCosScalar(a)     PetscCosReal(a)</font></strong>
<a name="line548">548: </a><strong><font color="#228B22">  #define PetscTanScalar(a)     PetscTanReal(a)</font></strong>
<a name="line549">549: </a><strong><font color="#228B22">  #define PetscAsinScalar(a)    PetscAsinReal(a)</font></strong>
<a name="line550">550: </a><strong><font color="#228B22">  #define PetscAcosScalar(a)    PetscAcosReal(a)</font></strong>
<a name="line551">551: </a><strong><font color="#228B22">  #define PetscAtanScalar(a)    PetscAtanReal(a)</font></strong>
<a name="line552">552: </a><strong><font color="#228B22">  #define PetscSinhScalar(a)    PetscSinhReal(a)</font></strong>
<a name="line553">553: </a><strong><font color="#228B22">  #define PetscCoshScalar(a)    PetscCoshReal(a)</font></strong>
<a name="line554">554: </a><strong><font color="#228B22">  #define PetscTanhScalar(a)    PetscTanhReal(a)</font></strong>
<a name="line555">555: </a><strong><font color="#228B22">  #define PetscAsinhScalar(a)   PetscAsinhReal(a)</font></strong>
<a name="line556">556: </a><strong><font color="#228B22">  #define PetscAcoshScalar(a)   PetscAcoshReal(a)</font></strong>
<a name="line557">557: </a><strong><font color="#228B22">  #define PetscAtanhScalar(a)   PetscAtanhReal(a)</font></strong>

<a name="line559">559: </a><font color="#A020F0">#endif </font><font color="#B22222">/* PETSC_USE_COMPLEX */</font><font color="#A020F0"></font>

<a name="line561">561: </a><font color="#B22222">/*</font>
<a name="line562">562: </a><font color="#B22222">   Certain objects may be created using either single or double precision.</font>
<a name="line563">563: </a><font color="#B22222">   This is currently not used.</font>
<a name="line564">564: </a><font color="#B22222">*/</font>
<a name="line565">565: </a><font color="#4169E1">typedef</font> <font color="#4169E1">enum</font> {
<a name="line566">566: </a>  PETSC_SCALAR_DOUBLE,
<a name="line567">567: </a>  PETSC_SCALAR_SINGLE,
<a name="line568">568: </a>  PETSC_SCALAR_LONG_DOUBLE,
<a name="line569">569: </a>  PETSC_SCALAR_HALF
<a name="line570">570: </a>} PetscScalarPrecision;

<a name="line572">572: </a><font color="#B22222">/*MC</font>
<a name="line573">573: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscAbs.html">PetscAbs</a> - Returns the absolute value of a number</font>

<a name="line575">575: </a><font color="#B22222">   Synopsis:</font>
<a name="line576">576: </a>#include <A href="../include/petscmath.h.html">&lt;petscmath.h&gt;</A>
<a name="line577">577: </a><font color="#B22222">   type <a href="../manualpages/Sys/PetscAbs.html">PetscAbs</a>(type v)</font>

<a name="line579">579: </a><font color="#B22222">   Not Collective</font>

<a name="line581">581: </a><font color="#B22222">   Input Parameter:</font>
<a name="line582">582: </a><font color="#B22222">.  v - the number</font>

<a name="line584">584: </a><font color="#B22222">   Level: beginner</font>

<a name="line586">586: </a><font color="#B22222">   Note:</font>
<a name="line587">587: </a><font color="#B22222">   The type can be integer or real floating point value, but cannot be complex</font>

<a name="line589">589: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscAbsInt.html">PetscAbsInt</a>()`, `<a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>()`, `PetscAbsScalar()`, `<a href="../manualpages/Sys/PetscSign.html">PetscSign</a>()`</font>
<a name="line590">590: </a><font color="#B22222">M*/</font>
<a name="line591">591: </a><strong><font color="#228B22">#define <a href="../manualpages/Sys/PetscAbs.html">PetscAbs</a>(a) (((a) &gt;= 0) ? (a) : (-(a)))</font></strong>

<a name="line593">593: </a><font color="#B22222">/*MC</font>
<a name="line594">594: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscSign.html">PetscSign</a> - Returns the sign of a number as an integer of value -1, 0, or 1</font>

<a name="line596">596: </a><font color="#B22222">   Synopsis:</font>
<a name="line597">597: </a>#include <A href="../include/petscmath.h.html">&lt;petscmath.h&gt;</A>
<a name="line598">598: </a><font color="#B22222">   int <a href="../manualpages/Sys/PetscSign.html">PetscSign</a>(type v)</font>

<a name="line600">600: </a><font color="#B22222">   Not Collective</font>

<a name="line602">602: </a><font color="#B22222">   Input Parameter:</font>
<a name="line603">603: </a><font color="#B22222">.  v - the number</font>

<a name="line605">605: </a><font color="#B22222">   Level: beginner</font>

<a name="line607">607: </a><font color="#B22222">   Note:</font>
<a name="line608">608: </a><font color="#B22222">   The type can be integer or real floating point value</font>

<a name="line610">610: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscAbsInt.html">PetscAbsInt</a>()`, `<a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>()`, `PetscAbsScalar()`</font>
<a name="line611">611: </a><font color="#B22222">M*/</font>
<a name="line612">612: </a><strong><font color="#228B22">#define <a href="../manualpages/Sys/PetscSign.html">PetscSign</a>(a) (((a) &gt;= 0) ? ((a) == 0 ? 0 : 1) : -1)</font></strong>

<a name="line614">614: </a><font color="#B22222">/*MC</font>
<a name="line615">615: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscMin.html">PetscMin</a> - Returns minimum of two numbers</font>

<a name="line617">617: </a><font color="#B22222">   Synopsis:</font>
<a name="line618">618: </a>#include <A href="../include/petscmath.h.html">&lt;petscmath.h&gt;</A>
<a name="line619">619: </a><font color="#B22222">   type <a href="../manualpages/Sys/PetscMin.html">PetscMin</a>(type v1,type v2)</font>

<a name="line621">621: </a><font color="#B22222">   Not Collective</font>

<a name="line623">623: </a><font color="#B22222">   Input Parameters:</font>
<a name="line624">624: </a><font color="#B22222">+  v1 - first value to find minimum of</font>
<a name="line625">625: </a><font color="#B22222">-  v2 - second value to find minimum of</font>

<a name="line627">627: </a><font color="#B22222">   Level: beginner</font>

<a name="line629">629: </a><font color="#B22222">   Note:</font>
<a name="line630">630: </a><font color="#B22222">   The type can be integer or floating point value, but cannot be complex</font>

<a name="line632">632: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscMax.html">PetscMax</a>()`, `<a href="../manualpages/Sys/PetscClipInterval.html">PetscClipInterval</a>()`, `<a href="../manualpages/Sys/PetscAbsInt.html">PetscAbsInt</a>()`, `<a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>()`, `<a href="../manualpages/Sys/PetscSqr.html">PetscSqr</a>()`</font>
<a name="line633">633: </a><font color="#B22222">M*/</font>
<a name="line634">634: </a><strong><font color="#228B22">#define <a href="../manualpages/Sys/PetscMin.html">PetscMin</a>(a, b) (((a) &lt; (b)) ? (a) : (b))</font></strong>

<a name="line636">636: </a><font color="#B22222">/*MC</font>
<a name="line637">637: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscMax.html">PetscMax</a> - Returns maximum of two numbers</font>

<a name="line639">639: </a><font color="#B22222">   Synopsis:</font>
<a name="line640">640: </a>#include <A href="../include/petscmath.h.html">&lt;petscmath.h&gt;</A>
<a name="line641">641: </a><font color="#B22222">   type max <a href="../manualpages/Sys/PetscMax.html">PetscMax</a>(type v1,type v2)</font>

<a name="line643">643: </a><font color="#B22222">   Not Collective</font>

<a name="line645">645: </a><font color="#B22222">   Input Parameters:</font>
<a name="line646">646: </a><font color="#B22222">+  v1 - first value to find maximum of</font>
<a name="line647">647: </a><font color="#B22222">-  v2 - second value to find maximum of</font>

<a name="line649">649: </a><font color="#B22222">   Level: beginner</font>

<a name="line651">651: </a><font color="#B22222">   Note:</font>
<a name="line652">652: </a><font color="#B22222">   The type can be integer or floating point value</font>

<a name="line654">654: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscMin.html">PetscMin</a>()`, `<a href="../manualpages/Sys/PetscClipInterval.html">PetscClipInterval</a>()`, `<a href="../manualpages/Sys/PetscAbsInt.html">PetscAbsInt</a>()`, `<a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>()`, `<a href="../manualpages/Sys/PetscSqr.html">PetscSqr</a>()`</font>
<a name="line655">655: </a><font color="#B22222">M*/</font>
<a name="line656">656: </a><strong><font color="#228B22">#define <a href="../manualpages/Sys/PetscMax.html">PetscMax</a>(a, b) (((a) &lt; (b)) ? (b) : (a))</font></strong>

<a name="line658">658: </a><font color="#B22222">/*MC</font>
<a name="line659">659: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscClipInterval.html">PetscClipInterval</a> - Returns a number clipped to be within an interval</font>

<a name="line661">661: </a><font color="#B22222">   Synopsis:</font>
<a name="line662">662: </a>#include <A href="../include/petscmath.h.html">&lt;petscmath.h&gt;</A>
<a name="line663">663: </a><font color="#B22222">   type clip <a href="../manualpages/Sys/PetscClipInterval.html">PetscClipInterval</a>(type x,type a,type b)</font>

<a name="line665">665: </a><font color="#B22222">   Not Collective</font>

<a name="line667">667: </a><font color="#B22222">   Input Parameters:</font>
<a name="line668">668: </a><font color="#B22222">+  x - value to use if within interval [a,b]</font>
<a name="line669">669: </a><font color="#B22222">.  a - lower end of interval</font>
<a name="line670">670: </a><font color="#B22222">-  b - upper end of interval</font>

<a name="line672">672: </a><font color="#B22222">   Level: beginner</font>

<a name="line674">674: </a><font color="#B22222">   Note:</font>
<a name="line675">675: </a><font color="#B22222">   The type can be integer or floating point value</font>

<a name="line677">677: </a><font color="#B22222">   Example\:</font>
<a name="line678">678: </a><font color="#B22222">.vb</font>
<a name="line679">679: </a><font color="#B22222">  <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> c = <a href="../manualpages/Sys/PetscClipInterval.html">PetscClipInterval</a>(5, 2, 3); // the value of c is 3</font>
<a name="line680">680: </a><font color="#B22222">  <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> c = <a href="../manualpages/Sys/PetscClipInterval.html">PetscClipInterval</a>(5, 2, 6); // the value of c is 5</font>
<a name="line681">681: </a><font color="#B22222">.ve</font>

<a name="line683">683: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscMin.html">PetscMin</a>()`, `<a href="../manualpages/Sys/PetscMax.html">PetscMax</a>()`, `<a href="../manualpages/Sys/PetscAbsInt.html">PetscAbsInt</a>()`, `<a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>()`, `<a href="../manualpages/Sys/PetscSqr.html">PetscSqr</a>()`</font>
<a name="line684">684: </a><font color="#B22222">M*/</font>
<a name="line685">685: </a><strong><font color="#228B22">#define <a href="../manualpages/Sys/PetscClipInterval.html">PetscClipInterval</a>(x, a, b) (<a href="../manualpages/Sys/PetscMax.html">PetscMax</a>((a), <a href="../manualpages/Sys/PetscMin.html">PetscMin</a>((x), (b))))</font></strong>

<a name="line687">687: </a><font color="#B22222">/*MC</font>
<a name="line688">688: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscAbsInt.html">PetscAbsInt</a> - Returns the absolute value of an integer</font>

<a name="line690">690: </a><font color="#B22222">   Synopsis:</font>
<a name="line691">691: </a>#include <A href="../include/petscmath.h.html">&lt;petscmath.h&gt;</A>
<a name="line692">692: </a><font color="#B22222">   int abs <a href="../manualpages/Sys/PetscAbsInt.html">PetscAbsInt</a>(int v1)</font>

<a name="line694">694: </a><font color="#B22222">   Input Parameter:</font>
<a name="line695">695: </a><font color="#B22222">.   v1 - the integer</font>

<a name="line697">697: </a><font color="#B22222">   Level: beginner</font>

<a name="line699">699: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscMax.html">PetscMax</a>()`, `<a href="../manualpages/Sys/PetscMin.html">PetscMin</a>()`, `<a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>()`, `<a href="../manualpages/Sys/PetscSqr.html">PetscSqr</a>()`</font>
<a name="line700">700: </a><font color="#B22222">M*/</font>
<a name="line701">701: </a><strong><font color="#228B22">#define <a href="../manualpages/Sys/PetscAbsInt.html">PetscAbsInt</a>(a) (((a) &lt; 0) ? (-(a)) : (a))</font></strong>

<a name="line703">703: </a><font color="#B22222">/*MC</font>
<a name="line704">704: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a> - Returns the absolute value of a real number</font>

<a name="line706">706: </a><font color="#B22222">   Synopsis:</font>
<a name="line707">707: </a>#include <A href="../include/petscmath.h.html">&lt;petscmath.h&gt;</A>
<a name="line708">708: </a><font color="#B22222">   Real abs <a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>(<a href="../manualpages/Sys/PetscReal.html">PetscReal</a> v1)</font>

<a name="line710">710: </a><font color="#B22222">   Input Parameter:</font>
<a name="line711">711: </a><font color="#B22222">.   v1 - the `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>` value</font>

<a name="line713">713: </a><font color="#B22222">   Level: beginner</font>

<a name="line715">715: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>`, `<a href="../manualpages/Sys/PetscMax.html">PetscMax</a>()`, `<a href="../manualpages/Sys/PetscMin.html">PetscMin</a>()`, `<a href="../manualpages/Sys/PetscAbsInt.html">PetscAbsInt</a>()`, `<a href="../manualpages/Sys/PetscSqr.html">PetscSqr</a>()`</font>
<a name="line716">716: </a><font color="#B22222">M*/</font>
<a name="line717">717: </a><font color="#A020F0">#if defined(PETSC_USE_REAL_SINGLE)</font>
<a name="line718">718: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>(a) fabsf(a)</font></strong>
<a name="line719">719: </a><font color="#A020F0">#elif defined(PETSC_USE_REAL_DOUBLE)</font>
<a name="line720">720: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>(a) fabs(a)</font></strong>
<a name="line721">721: </a><font color="#A020F0">#elif defined(PETSC_USE_REAL___FLOAT128)</font>
<a name="line722">722: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>(a) fabsq(a)</font></strong>
<a name="line723">723: </a><font color="#A020F0">#elif defined(PETSC_USE_REAL___FP16)</font>
<a name="line724">724: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>(a) fabsf(a)</font></strong>
<a name="line725">725: </a><font color="#A020F0">#endif</font>

<a name="line727">727: </a><font color="#B22222">/*MC</font>
<a name="line728">728: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscSqr.html">PetscSqr</a> - Returns the square of a number</font>

<a name="line730">730: </a><font color="#B22222">   Synopsis:</font>
<a name="line731">731: </a>#include <A href="../include/petscmath.h.html">&lt;petscmath.h&gt;</A>
<a name="line732">732: </a><font color="#B22222">   type sqr <a href="../manualpages/Sys/PetscSqr.html">PetscSqr</a>(type v1)</font>

<a name="line734">734: </a><font color="#B22222">   Not Collective</font>

<a name="line736">736: </a><font color="#B22222">   Input Parameter:</font>
<a name="line737">737: </a><font color="#B22222">.   v1 - the value</font>

<a name="line739">739: </a><font color="#B22222">   Level: beginner</font>

<a name="line741">741: </a><font color="#B22222">   Note:</font>
<a name="line742">742: </a><font color="#B22222">   The type can be integer, floating point, or complex floating point</font>

<a name="line744">744: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscMax.html">PetscMax</a>()`, `<a href="../manualpages/Sys/PetscMin.html">PetscMin</a>()`, `<a href="../manualpages/Sys/PetscAbsInt.html">PetscAbsInt</a>()`, `<a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>()`</font>
<a name="line745">745: </a><font color="#B22222">M*/</font>
<a name="line746">746: </a><strong><font color="#228B22">#define <a href="../manualpages/Sys/PetscSqr.html">PetscSqr</a>(a) ((a) * (a))</font></strong>

<a name="line748">748: </a><font color="#B22222">/*MC</font>
<a name="line749">749: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscRealConstant.html">PetscRealConstant</a> - a compile time macro that ensures a given constant real number is properly represented in the configured</font>
<a name="line750">750: </a><font color="#B22222">   precision of `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>` be it half, single, double or 128-bit representation</font>

<a name="line752">752: </a><font color="#B22222">   Synopsis:</font>
<a name="line753">753: </a>#include <A href="../include/petscmath.h.html">&lt;petscmath.h&gt;</A>
<a name="line754">754: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> <a href="../manualpages/Sys/PetscRealConstant.html">PetscRealConstant</a>(real_number)</font>

<a name="line756">756: </a><font color="#B22222">   Not Collective</font>

<a name="line758">758: </a><font color="#B22222">   Input Parameter:</font>
<a name="line759">759: </a><font color="#B22222">.   v1 - the real number, for example 1.5</font>

<a name="line761">761: </a><font color="#B22222">   Level: beginner</font>

<a name="line763">763: </a><font color="#B22222">   Note:</font>
<a name="line764">764: </a><font color="#B22222">   For example, if PETSc is configured with `--with-precision=__float128` and one writes</font>
<a name="line765">765: </a><font color="#B22222">.vb</font>
<a name="line766">766: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> d = 1.5;</font>
<a name="line767">767: </a><font color="#B22222">.ve</font>
<a name="line768">768: </a><font color="#B22222">   the result is 1.5 in double precision extended to 128 bit representation, meaning it is very far from the correct value. Hence, one should write</font>
<a name="line769">769: </a><font color="#B22222">.vb</font>
<a name="line770">770: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> d = <a href="../manualpages/Sys/PetscRealConstant.html">PetscRealConstant</a>(1.5);</font>
<a name="line771">771: </a><font color="#B22222">.ve</font>

<a name="line773">773: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>`</font>
<a name="line774">774: </a><font color="#B22222">M*/</font>
<a name="line775">775: </a><font color="#A020F0">#if defined(PETSC_USE_REAL_SINGLE)</font>
<a name="line776">776: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PetscRealConstant.html">PetscRealConstant</a>(constant) constant##F</font></strong>
<a name="line777">777: </a><font color="#A020F0">#elif defined(PETSC_USE_REAL_DOUBLE)</font>
<a name="line778">778: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PetscRealConstant.html">PetscRealConstant</a>(constant) constant</font></strong>
<a name="line779">779: </a><font color="#A020F0">#elif defined(PETSC_USE_REAL___FLOAT128)</font>
<a name="line780">780: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PetscRealConstant.html">PetscRealConstant</a>(constant) constant##Q</font></strong>
<a name="line781">781: </a><font color="#A020F0">#elif defined(PETSC_USE_REAL___FP16)</font>
<a name="line782">782: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PetscRealConstant.html">PetscRealConstant</a>(constant) constant##F</font></strong>
<a name="line783">783: </a><font color="#A020F0">#endif</font>

<a name="line785">785: </a><font color="#B22222">/*</font>
<a name="line786">786: </a><font color="#B22222">     Basic constants</font>
<a name="line787">787: </a><font color="#B22222">*/</font>
<a name="line788">788: </a><font color="#B22222">/*MC</font>
<a name="line789">789: </a><font color="#B22222">  <a href="../manualpages/Sys/PETSC_PI.html">PETSC_PI</a> - the value of $ \pi$ to the correct precision of `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>`.</font>

<a name="line791">791: </a><font color="#B22222">  Level: beginner</font>

<a name="line793">793: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>`, `<a href="../manualpages/Sys/PETSC_PHI.html">PETSC_PHI</a>`, `<a href="../manualpages/Sys/PETSC_SQRT2.html">PETSC_SQRT2</a>`</font>
<a name="line794">794: </a><font color="#B22222">M*/</font>

<a name="line796">796: </a><font color="#B22222">/*MC</font>
<a name="line797">797: </a><font color="#B22222">  <a href="../manualpages/Sys/PETSC_PHI.html">PETSC_PHI</a> - the value of $ \phi$, the Golden Ratio, to the correct precision of `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>`.</font>

<a name="line799">799: </a><font color="#B22222">  Level: beginner</font>

<a name="line801">801: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>`, `<a href="../manualpages/Sys/PETSC_PI.html">PETSC_PI</a>`, `<a href="../manualpages/Sys/PETSC_SQRT2.html">PETSC_SQRT2</a>`</font>
<a name="line802">802: </a><font color="#B22222">M*/</font>

<a name="line804">804: </a><font color="#B22222">/*MC</font>
<a name="line805">805: </a><font color="#B22222">  <a href="../manualpages/Sys/PETSC_SQRT2.html">PETSC_SQRT2</a> - the value of $ \sqrt{2} $ to the correct precision of `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>`.</font>

<a name="line807">807: </a><font color="#B22222">  Level: beginner</font>

<a name="line809">809: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>`, `<a href="../manualpages/Sys/PETSC_PI.html">PETSC_PI</a>`, `<a href="../manualpages/Sys/PETSC_PHI.html">PETSC_PHI</a>`</font>
<a name="line810">810: </a><font color="#B22222">M*/</font>

<a name="line812">812: </a><strong><font color="#228B22">#define <a href="../manualpages/Sys/PETSC_PI.html">PETSC_PI</a>    <a href="../manualpages/Sys/PetscRealConstant.html">PetscRealConstant</a>(3.1415926535897932384626433832795029)</font></strong>
<a name="line813">813: </a><strong><font color="#228B22">#define <a href="../manualpages/Sys/PETSC_PHI.html">PETSC_PHI</a>   <a href="../manualpages/Sys/PetscRealConstant.html">PetscRealConstant</a>(1.6180339887498948482045868343656381)</font></strong>
<a name="line814">814: </a><strong><font color="#228B22">#define <a href="../manualpages/Sys/PETSC_SQRT2.html">PETSC_SQRT2</a> <a href="../manualpages/Sys/PetscRealConstant.html">PetscRealConstant</a>(1.4142135623730950488016887242096981)</font></strong>

<a name="line816">816: </a><font color="#B22222">/*MC</font>
<a name="line817">817: </a><font color="#B22222">  <a href="../manualpages/Sys/PETSC_MAX_REAL.html">PETSC_MAX_REAL</a> - the largest real value that can be stored in a `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>`</font>

<a name="line819">819: </a><font color="#B22222">  Level: beginner</font>

<a name="line821">821: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PETSC_MIN_REAL.html">PETSC_MIN_REAL</a>`, `<a href="../manualpages/Sys/PETSC_REAL_MIN.html">PETSC_REAL_MIN</a>`, `<a href="../manualpages/Sys/PETSC_MACHINE_EPSILON.html">PETSC_MACHINE_EPSILON</a>`, `<a href="../manualpages/Sys/PETSC_SQRT_MACHINE_EPSILON.html">PETSC_SQRT_MACHINE_EPSILON</a>`, `<a href="../manualpages/Sys/PETSC_SMALL.html">PETSC_SMALL</a>`</font>
<a name="line822">822: </a><font color="#B22222">M*/</font>

<a name="line824">824: </a><font color="#B22222">/*MC</font>
<a name="line825">825: </a><font color="#B22222">  <a href="../manualpages/Sys/PETSC_MIN_REAL.html">PETSC_MIN_REAL</a> - the smallest real value that can be stored in a `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>`, generally this is - `<a href="../manualpages/Sys/PETSC_MAX_REAL.html">PETSC_MAX_REAL</a>`</font>

<a name="line827">827: </a><font color="#B22222">  Level: beginner</font>

<a name="line829">829: </a><font color="#B22222">.seealso `<a href="../manualpages/Sys/PETSC_MAX_REAL.html">PETSC_MAX_REAL</a>`, `<a href="../manualpages/Sys/PETSC_REAL_MIN.html">PETSC_REAL_MIN</a>`, `<a href="../manualpages/Sys/PETSC_MACHINE_EPSILON.html">PETSC_MACHINE_EPSILON</a>`, `<a href="../manualpages/Sys/PETSC_SQRT_MACHINE_EPSILON.html">PETSC_SQRT_MACHINE_EPSILON</a>`, `<a href="../manualpages/Sys/PETSC_SMALL.html">PETSC_SMALL</a>`</font>
<a name="line830">830: </a><font color="#B22222">M*/</font>

<a name="line832">832: </a><font color="#B22222">/*MC</font>
<a name="line833">833: </a><font color="#B22222">  <a href="../manualpages/Sys/PETSC_REAL_MIN.html">PETSC_REAL_MIN</a> - the smallest positive normalized real value that can be stored in a `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>`.</font>

<a name="line835">835: </a><font color="#B22222">  Level: beginner</font>

<a name="line837">837: </a><font color="#B22222">  Note:</font>
<a name="line838">838: </a><font color="#B22222">  See &lt;https://en.wikipedia.org/wiki/Subnormal_number&gt; for a discussion of normalized and subnormal floating point numbers</font>

<a name="line840">840: </a><font color="#B22222">  Developer Note:</font>
<a name="line841">841: </a><font color="#B22222">  The naming is confusing as there is both a `<a href="../manualpages/Sys/PETSC_REAL_MIN.html">PETSC_REAL_MIN</a>` and `<a href="../manualpages/Sys/PETSC_MIN_REAL.html">PETSC_MIN_REAL</a>` with different meanings.</font>

<a name="line843">843: </a><font color="#B22222">.seealso `<a href="../manualpages/Sys/PETSC_MAX_REAL.html">PETSC_MAX_REAL</a>`, `<a href="../manualpages/Sys/PETSC_MIN_REAL.html">PETSC_MIN_REAL</a>`, `<a href="../manualpages/Sys/PETSC_MACHINE_EPSILON.html">PETSC_MACHINE_EPSILON</a>`, `<a href="../manualpages/Sys/PETSC_SQRT_MACHINE_EPSILON.html">PETSC_SQRT_MACHINE_EPSILON</a>`, `<a href="../manualpages/Sys/PETSC_SMALL.html">PETSC_SMALL</a>`</font>
<a name="line844">844: </a><font color="#B22222">M*/</font>

<a name="line846">846: </a><font color="#B22222">/*MC</font>
<a name="line847">847: </a><font color="#B22222">  <a href="../manualpages/Sys/PETSC_MACHINE_EPSILON.html">PETSC_MACHINE_EPSILON</a> - the machine epsilon for the precision of `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>`</font>

<a name="line849">849: </a><font color="#B22222">  Level: beginner</font>

<a name="line851">851: </a><font color="#B22222">  Note:</font>
<a name="line852">852: </a><font color="#B22222">  See &lt;https://en.wikipedia.org/wiki/Machine_epsilon&gt;</font>

<a name="line854">854: </a><font color="#B22222">.seealso `<a href="../manualpages/Sys/PETSC_MAX_REAL.html">PETSC_MAX_REAL</a>`, `<a href="../manualpages/Sys/PETSC_MIN_REAL.html">PETSC_MIN_REAL</a>`, `<a href="../manualpages/Sys/PETSC_REAL_MIN.html">PETSC_REAL_MIN</a>`, `<a href="../manualpages/Sys/PETSC_SQRT_MACHINE_EPSILON.html">PETSC_SQRT_MACHINE_EPSILON</a>`, `<a href="../manualpages/Sys/PETSC_SMALL.html">PETSC_SMALL</a>`</font>
<a name="line855">855: </a><font color="#B22222">M*/</font>

<a name="line857">857: </a><font color="#B22222">/*MC</font>
<a name="line858">858: </a><font color="#B22222">  <a href="../manualpages/Sys/PETSC_SQRT_MACHINE_EPSILON.html">PETSC_SQRT_MACHINE_EPSILON</a> - the square root of the machine epsilon for the precision of `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>`</font>

<a name="line860">860: </a><font color="#B22222">  Level: beginner</font>

<a name="line862">862: </a><font color="#B22222">  Note:</font>
<a name="line863">863: </a><font color="#B22222">  See `<a href="../manualpages/Sys/PETSC_MACHINE_EPSILON.html">PETSC_MACHINE_EPSILON</a>`</font>

<a name="line865">865: </a><font color="#B22222">.seealso `<a href="../manualpages/Sys/PETSC_MAX_REAL.html">PETSC_MAX_REAL</a>`, `<a href="../manualpages/Sys/PETSC_MIN_REAL.html">PETSC_MIN_REAL</a>`, `<a href="../manualpages/Sys/PETSC_REAL_MIN.html">PETSC_REAL_MIN</a>`, `<a href="../manualpages/Sys/PETSC_MACHINE_EPSILON.html">PETSC_MACHINE_EPSILON</a>`, `<a href="../manualpages/Sys/PETSC_SMALL.html">PETSC_SMALL</a>`</font>
<a name="line866">866: </a><font color="#B22222">M*/</font>

<a name="line868">868: </a><font color="#B22222">/*MC</font>
<a name="line869">869: </a><font color="#B22222">  <a href="../manualpages/Sys/PETSC_SMALL.html">PETSC_SMALL</a> - an arbitrary "small" number which depends on the precision of `<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>` used in some PETSc examples</font>
<a name="line870">870: </a><font color="#B22222">  and in `<a href="../manualpages/Sys/PetscApproximateLTE.html">PetscApproximateLTE</a>()` and `<a href="../manualpages/Sys/PetscApproximateGTE.html">PetscApproximateGTE</a>()` to determine if a computation was successful.</font>

<a name="line872">872: </a><font color="#B22222">  Level: beginner</font>

<a name="line874">874: </a><font color="#B22222">  Note:</font>
<a name="line875">875: </a><font color="#B22222">  See `<a href="../manualpages/Sys/PETSC_MACHINE_EPSILON.html">PETSC_MACHINE_EPSILON</a>`</font>

<a name="line877">877: </a><font color="#B22222">.seealso `<a href="../manualpages/Sys/PetscApproximateLTE.html">PetscApproximateLTE</a>()`, `<a href="../manualpages/Sys/PetscApproximateGTE.html">PetscApproximateGTE</a>()`, `<a href="../manualpages/Sys/PETSC_MAX_REAL.html">PETSC_MAX_REAL</a>`, `<a href="../manualpages/Sys/PETSC_MIN_REAL.html">PETSC_MIN_REAL</a>`, `<a href="../manualpages/Sys/PETSC_REAL_MIN.html">PETSC_REAL_MIN</a>`, `<a href="../manualpages/Sys/PETSC_MACHINE_EPSILON.html">PETSC_MACHINE_EPSILON</a>`,</font>
<a name="line878">878: </a><font color="#B22222">         `<a href="../manualpages/Sys/PETSC_SQRT_MACHINE_EPSILON.html">PETSC_SQRT_MACHINE_EPSILON</a>`</font>
<a name="line879">879: </a><font color="#B22222">M*/</font>

<a name="line881">881: </a><font color="#A020F0">#if defined(PETSC_USE_REAL_SINGLE)</font>
<a name="line882">882: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_MAX_REAL.html">PETSC_MAX_REAL</a>             3.40282346638528860e+38F</font></strong>
<a name="line883">883: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_MIN_REAL.html">PETSC_MIN_REAL</a>             (-<a href="../manualpages/Sys/PETSC_MAX_REAL.html">PETSC_MAX_REAL</a>)</font></strong>
<a name="line884">884: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_REAL_MIN.html">PETSC_REAL_MIN</a>             1.1754944e-38F</font></strong>
<a name="line885">885: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_MACHINE_EPSILON.html">PETSC_MACHINE_EPSILON</a>      1.19209290e-07F</font></strong>
<a name="line886">886: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_SQRT_MACHINE_EPSILON.html">PETSC_SQRT_MACHINE_EPSILON</a> 3.45266983e-04F</font></strong>
<a name="line887">887: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_SMALL.html">PETSC_SMALL</a>                1.e-5F</font></strong>
<a name="line888">888: </a><font color="#A020F0">#elif defined(PETSC_USE_REAL_DOUBLE)</font>
<a name="line889">889: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_MAX_REAL.html">PETSC_MAX_REAL</a>             1.7976931348623157e+308</font></strong>
<a name="line890">890: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_MIN_REAL.html">PETSC_MIN_REAL</a>             (-<a href="../manualpages/Sys/PETSC_MAX_REAL.html">PETSC_MAX_REAL</a>)</font></strong>
<a name="line891">891: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_REAL_MIN.html">PETSC_REAL_MIN</a>             2.225073858507201e-308</font></strong>
<a name="line892">892: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_MACHINE_EPSILON.html">PETSC_MACHINE_EPSILON</a>      2.2204460492503131e-16</font></strong>
<a name="line893">893: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_SQRT_MACHINE_EPSILON.html">PETSC_SQRT_MACHINE_EPSILON</a> 1.490116119384766e-08</font></strong>
<a name="line894">894: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_SMALL.html">PETSC_SMALL</a>                1.e-10</font></strong>
<a name="line895">895: </a><font color="#A020F0">#elif defined(PETSC_USE_REAL___FLOAT128)</font>
<a name="line896">896: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_MAX_REAL.html">PETSC_MAX_REAL</a>             FLT128_MAX</font></strong>
<a name="line897">897: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_MIN_REAL.html">PETSC_MIN_REAL</a>             (-FLT128_MAX)</font></strong>
<a name="line898">898: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_REAL_MIN.html">PETSC_REAL_MIN</a>             FLT128_MIN</font></strong>
<a name="line899">899: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_MACHINE_EPSILON.html">PETSC_MACHINE_EPSILON</a>      FLT128_EPSILON</font></strong>
<a name="line900">900: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_SQRT_MACHINE_EPSILON.html">PETSC_SQRT_MACHINE_EPSILON</a> 1.38777878078144567552953958511352539e-17Q</font></strong>
<a name="line901">901: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_SMALL.html">PETSC_SMALL</a>                1.e-20Q</font></strong>
<a name="line902">902: </a><font color="#A020F0">#elif defined(PETSC_USE_REAL___FP16)</font>
<a name="line903">903: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_MAX_REAL.html">PETSC_MAX_REAL</a>             65504.0F</font></strong>
<a name="line904">904: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_MIN_REAL.html">PETSC_MIN_REAL</a>             (-<a href="../manualpages/Sys/PETSC_MAX_REAL.html">PETSC_MAX_REAL</a>)</font></strong>
<a name="line905">905: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_REAL_MIN.html">PETSC_REAL_MIN</a>             .00006103515625F</font></strong>
<a name="line906">906: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_MACHINE_EPSILON.html">PETSC_MACHINE_EPSILON</a>      .0009765625F</font></strong>
<a name="line907">907: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_SQRT_MACHINE_EPSILON.html">PETSC_SQRT_MACHINE_EPSILON</a> .03125F</font></strong>
<a name="line908">908: </a><strong><font color="#228B22">  #define <a href="../manualpages/Sys/PETSC_SMALL.html">PETSC_SMALL</a>                5.e-3F</font></strong>
<a name="line909">909: </a><font color="#A020F0">#endif</font>

<a name="line911">911: </a><font color="#B22222">/*MC</font>
<a name="line912">912: </a><font color="#B22222">  <a href="../manualpages/Sys/PETSC_INFINITY.html">PETSC_INFINITY</a> - a finite number that represents infinity for setting certain bounds in `<a href="../manualpages/Tao/Tao.html">Tao</a>`</font>

<a name="line914">914: </a><font color="#B22222">  Level: intermediate</font>

<a name="line916">916: </a><font color="#B22222">  Note:</font>
<a name="line917">917: </a><font color="#B22222">  This is not the IEEE infinity value</font>

<a name="line919">919: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PETSC_NINFINITY.html">PETSC_NINFINITY</a>`, `<a href="../manualpages/SNES/SNESVIGetVariableBounds.html">SNESVIGetVariableBounds</a>()`, `<a href="../manualpages/SNES/SNESVISetComputeVariableBounds.html">SNESVISetComputeVariableBounds</a>()`, `<a href="../manualpages/SNES/SNESVISetVariableBounds.html">SNESVISetVariableBounds</a>()`</font>
<a name="line920">920: </a><font color="#B22222">M*/</font>
<a name="line921">921: </a><strong><font color="#228B22">#define <a href="../manualpages/Sys/PETSC_INFINITY.html">PETSC_INFINITY</a> (<a href="../manualpages/Sys/PETSC_MAX_REAL.html">PETSC_MAX_REAL</a> / 4)</font></strong>

<a name="line923">923: </a><font color="#B22222">/*MC</font>
<a name="line924">924: </a><font color="#B22222">  <a href="../manualpages/Sys/PETSC_NINFINITY.html">PETSC_NINFINITY</a> - a finite number that represents negative infinity for setting certain bounds in `<a href="../manualpages/Tao/Tao.html">Tao</a>`</font>

<a name="line926">926: </a><font color="#B22222">  Level: intermediate</font>

<a name="line928">928: </a><font color="#B22222">  Note:</font>
<a name="line929">929: </a><font color="#B22222">  This is not the negative IEEE infinity value</font>

<a name="line931">931: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PETSC_INFINITY.html">PETSC_INFINITY</a>`, `<a href="../manualpages/SNES/SNESVIGetVariableBounds.html">SNESVIGetVariableBounds</a>()`, `<a href="../manualpages/SNES/SNESVISetComputeVariableBounds.html">SNESVISetComputeVariableBounds</a>()`, `<a href="../manualpages/SNES/SNESVISetVariableBounds.html">SNESVISetVariableBounds</a>()`</font>
<a name="line932">932: </a><font color="#B22222">M*/</font>
<a name="line933">933: </a><strong><font color="#228B22">#define <a href="../manualpages/Sys/PETSC_NINFINITY.html">PETSC_NINFINITY</a> (-<a href="../manualpages/Sys/PETSC_INFINITY.html">PETSC_INFINITY</a>)</font></strong>

<a name="line935">935: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscBool.html">PetscBool</a>  <a href="../manualpages/Sys/PetscIsInfReal.html">PetscIsInfReal</a>(<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)</font></strong>;
<a name="line936">936: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscBool.html">PetscBool</a>  <a href="../manualpages/Sys/PetscIsNanReal.html">PetscIsNanReal</a>(<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)</font></strong>;
<a name="line937">937: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscBool.html">PetscBool</a>  <a href="../manualpages/Sys/PetscIsNormalReal.html">PetscIsNormalReal</a>(<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)</font></strong>;
<a name="line938">938: </a><strong><font color="#4169E1"><a name="PetscIsInfOrNanReal"></a>static inline <a href="../manualpages/Sys/PetscBool.html">PetscBool</a> PetscIsInfOrNanReal(<a href="../manualpages/Sys/PetscReal.html">PetscReal</a> v)</font></strong>
<a name="line939">939: </a>{
<a name="line940">940: </a>  <font color="#4169E1">return</font> <a href="../manualpages/Sys/PetscIsInfReal.html">PetscIsInfReal</a>(v) || <a href="../manualpages/Sys/PetscIsNanReal.html">PetscIsNanReal</a>(v) ? <a href="../manualpages/Sys/PETSC_TRUE.html">PETSC_TRUE</a> : <a href="../manualpages/Sys/PETSC_FALSE.html">PETSC_FALSE</a>;
<a name="line941">941: </a>}
<a name="line942">942: </a><strong><font color="#4169E1"><a name="PetscIsInfScalar"></a>static inline <a href="../manualpages/Sys/PetscBool.html">PetscBool</a> PetscIsInfScalar(<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> v)</font></strong>
<a name="line943">943: </a>{
<a name="line944">944: </a>  <font color="#4169E1">return</font> <a href="../manualpages/Sys/PetscIsInfReal.html">PetscIsInfReal</a>(PetscAbsScalar(v));
<a name="line945">945: </a>}
<a name="line946">946: </a><strong><font color="#4169E1"><a name="PetscIsNanScalar"></a>static inline <a href="../manualpages/Sys/PetscBool.html">PetscBool</a> PetscIsNanScalar(<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> v)</font></strong>
<a name="line947">947: </a>{
<a name="line948">948: </a>  <font color="#4169E1">return</font> <a href="../manualpages/Sys/PetscIsNanReal.html">PetscIsNanReal</a>(PetscAbsScalar(v));
<a name="line949">949: </a>}
<a name="line950">950: </a><strong><font color="#4169E1"><a name="PetscIsInfOrNanScalar"></a>static inline <a href="../manualpages/Sys/PetscBool.html">PetscBool</a> PetscIsInfOrNanScalar(<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> v)</font></strong>
<a name="line951">951: </a>{
<a name="line952">952: </a>  <font color="#4169E1">return</font> PetscIsInfOrNanReal(PetscAbsScalar(v));
<a name="line953">953: </a>}
<a name="line954">954: </a><strong><font color="#4169E1"><a name="PetscIsNormalScalar"></a>static inline <a href="../manualpages/Sys/PetscBool.html">PetscBool</a> PetscIsNormalScalar(<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> v)</font></strong>
<a name="line955">955: </a>{
<a name="line956">956: </a>  <font color="#4169E1">return</font> <a href="../manualpages/Sys/PetscIsNormalReal.html">PetscIsNormalReal</a>(PetscAbsScalar(v));
<a name="line957">957: </a>}

<a name="line959">959: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscBool.html">PetscBool</a> <a href="../manualpages/Sys/PetscIsCloseAtTol.html">PetscIsCloseAtTol</a>(<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)</font></strong>;
<a name="line960">960: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscBool.html">PetscBool</a> <a href="../manualpages/Sys/PetscEqualReal.html">PetscEqualReal</a>(<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)</font></strong>;
<a name="line961">961: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscBool.html">PetscBool</a> <a href="../manualpages/Sys/PetscEqualScalar.html">PetscEqualScalar</a>(<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a>, <a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a>)</font></strong>;

<a name="line963">963: </a><font color="#B22222">/*@C</font>
<a name="line964">964: </a><font color="#B22222">  <a href="../manualpages/Sys/PetscIsCloseAtTolScalar.html">PetscIsCloseAtTolScalar</a> - Like `<a href="../manualpages/Sys/PetscIsCloseAtTol.html">PetscIsCloseAtTol</a>()` but for `<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a>`</font>

<a name="line966">966: </a><font color="#B22222">  Input Parameters:</font>
<a name="line967">967: </a><font color="#B22222">+ lhs  - The first number</font>
<a name="line968">968: </a><font color="#B22222">. rhs  - The second number</font>
<a name="line969">969: </a><font color="#B22222">. rtol - The relative tolerance</font>
<a name="line970">970: </a><font color="#B22222">- atol - The absolute tolerance</font>

<a name="line972">972: </a><font color="#B22222">  Level: beginner</font>

<a name="line974">974: </a><font color="#B22222">  Note:</font>
<a name="line975">975: </a><font color="#B22222">  This routine is equivalent to `<a href="../manualpages/Sys/PetscIsCloseAtTol.html">PetscIsCloseAtTol</a>()` when PETSc is configured without complex</font>
<a name="line976">976: </a><font color="#B22222">  numbers.</font>

<a name="line978">978: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscIsCloseAtTol.html">PetscIsCloseAtTol</a>()`</font>
<a name="line979">979: </a><font color="#B22222">@*/</font>
<a name="line980">980: </a><strong><font color="#4169E1"><a name="PetscIsCloseAtTolScalar"></a>static inline <a href="../manualpages/Sys/PetscBool.html">PetscBool</a> <a href="../manualpages/Sys/PetscIsCloseAtTolScalar.html">PetscIsCloseAtTolScalar</a>(<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> lhs, <a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> rhs, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> rtol, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> atol)</font></strong>
<a name="line981">981: </a>{
<a name="line982">982: </a>  <a href="../manualpages/Sys/PetscBool.html">PetscBool</a> close = <a href="../manualpages/Sys/PetscIsCloseAtTol.html">PetscIsCloseAtTol</a>(<a href="../manualpages/Sys/PetscRealPart.html">PetscRealPart</a>(lhs), <a href="../manualpages/Sys/PetscRealPart.html">PetscRealPart</a>(rhs), rtol, atol);

<a name="line984">984: </a>  <font color="#4169E1">if</font> (<a href="../manualpages/Sys/PetscDefined.html">PetscDefined</a>(USE_COMPLEX)) close = (<a href="../manualpages/Sys/PetscBool.html">PetscBool</a>)(close &amp;&amp; <a href="../manualpages/Sys/PetscIsCloseAtTol.html">PetscIsCloseAtTol</a>(<a href="../manualpages/Sys/PetscImaginaryPart.html">PetscImaginaryPart</a>(lhs), <a href="../manualpages/Sys/PetscImaginaryPart.html">PetscImaginaryPart</a>(rhs), rtol, atol));
<a name="line985">985: </a>  <font color="#4169E1">return</font> close;
<a name="line986">986: </a>}

<a name="line988">988: </a><font color="#B22222">/*</font>
<a name="line989">989: </a><font color="#B22222">    These macros are currently hardwired to match the regular data types, so there is no support for a different</font>
<a name="line990">990: </a><font color="#B22222">    MatScalar from <a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a>. We left the MatScalar in the source just in case we use it again.</font>
<a name="line991">991: </a><font color="#B22222"> */</font>
<a name="line992">992: </a><strong><font color="#228B22">#define MPIU_MATSCALAR <a href="../manualpages/Sys/MPIU_SCALAR.html">MPIU_SCALAR</a></font></strong>
<a name="line993">993: </a><font color="#4169E1">typedef <a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> MatScalar;</font>
<a name="line994">994: </a><font color="#4169E1">typedef <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>   MatReal;</font>

<a name="line996">996: </a><font color="#4169E1"><a name="petsc_mpiu_2scalar"></a>struct petsc_mpiu_2scalar </font>{
<a name="line997">997: </a>  <a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> a, b;
<a name="line998">998: </a>};
<a name="line999">999: </a><strong><font color="#4169E1">PETSC_EXTERN MPI_Datatype MPIU_2SCALAR PETSC_ATTRIBUTE_MPI_TYPE_TAG_LAYOUT_COMPATIBLE(struct petsc_mpiu_2scalar)</font></strong>;

<a name="line1001">1001: </a><font color="#B22222">/* MPI Datatypes for composite reductions */</font>
<a name="line1002">1002: </a><font color="#4169E1"><a name="petsc_mpiu_real_int"></a>struct petsc_mpiu_real_int </font>{
<a name="line1003">1003: </a>  <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> v;
<a name="line1004">1004: </a>  <a href="../manualpages/Sys/PetscInt.html">PetscInt</a>  i;
<a name="line1005">1005: </a>};

<a name="line1007">1007: </a><font color="#4169E1"><a name="petsc_mpiu_scalar_int"></a>struct petsc_mpiu_scalar_int </font>{
<a name="line1008">1008: </a>  <a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> v;
<a name="line1009">1009: </a>  <a href="../manualpages/Sys/PetscInt.html">PetscInt</a>    i;
<a name="line1010">1010: </a>};

<a name="line1012">1012: </a><strong><font color="#4169E1">PETSC_EXTERN MPI_Datatype MPIU_REAL_INT PETSC_ATTRIBUTE_MPI_TYPE_TAG_LAYOUT_COMPATIBLE(struct petsc_mpiu_real_int)</font></strong>;
<a name="line1013">1013: </a><strong><font color="#4169E1">PETSC_EXTERN MPI_Datatype MPIU_SCALAR_INT PETSC_ATTRIBUTE_MPI_TYPE_TAG_LAYOUT_COMPATIBLE(struct petsc_mpiu_scalar_int)</font></strong>;

<a name="line1015">1015: </a><font color="#A020F0">#if defined(PETSC_USE_64BIT_INDICES)</font>
<a name="line1016">1016: </a><font color="#4169E1">struct</font> <font color="#B22222">/* __attribute__((packed, aligned(alignof(<a href="../manualpages/Sys/PetscInt.html">PetscInt</a> *)))) */</font> petsc_mpiu_2int {
<a name="line1017">1017: </a>  <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> a;
<a name="line1018">1018: </a>  <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> b;
<a name="line1019">1019: </a>};
<a name="line1020">1020: </a><strong><font color="#4169E1"><a name="__attribute__"></a>struct __attribute__((packed))</font></strong> petsc_mpiu_int_mpiint {
<a name="line1021">1021: </a>  <a href="../manualpages/Sys/PetscInt.html">PetscInt</a>    a;
<a name="line1022">1022: </a>  <a href="../manualpages/Sys/PetscMPIInt.html">PetscMPIInt</a> b;
<a name="line1023">1023: </a>};
<a name="line1024">1024: </a><font color="#B22222">/*</font>
<a name="line1025">1025: </a><font color="#B22222"> static_assert(sizeof(struct petsc_mpiu_2int) == 2 * sizeof(<a href="../manualpages/Sys/PetscInt.html">PetscInt</a>), "");</font>
<a name="line1026">1026: </a><font color="#B22222"> static_assert(alignof(struct petsc_mpiu_2int) == alignof(<a href="../manualpages/Sys/PetscInt.html">PetscInt</a> *), "");</font>
<a name="line1027">1027: </a><font color="#B22222"> static_assert(alignof(struct petsc_mpiu_2int) == alignof(<a href="../manualpages/Sys/PetscInt.html">PetscInt</a>[2]), "");</font>

<a name="line1029">1029: </a><font color="#B22222"> clang generates warnings that petsc_mpiu_2int is not layout compatible with <a href="../manualpages/Sys/PetscInt.html">PetscInt</a>[2] or</font>
<a name="line1030">1030: </a><font color="#B22222"> <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> *, even though (with everything else uncommented) both of the static_asserts above</font>
<a name="line1031">1031: </a><font color="#B22222"> pass! So we just comment it out...</font>
<a name="line1032">1032: </a><font color="#B22222">*/</font>
<a name="line1033">1033: </a>PETSC_EXTERN MPI_Datatype MPIU_2INT <font color="#B22222">/* PETSC_ATTRIBUTE_MPI_TYPE_TAG_LAYOUT_COMPATIBLE(struct petsc_mpiu_2int) */</font>;
<a name="line1034">1034: </a>PETSC_EXTERN MPI_Datatype MPIU_INT_MPIINT <font color="#B22222">/* PETSC_ATTRIBUTE_MPI_TYPE_TAG_LAYOUT_COMPATIBLE(struct petsc_mpiu_int_mpiint) */</font>;
<a name="line1035">1035: </a><font color="#A020F0">#else</font>
<a name="line1036">1036: </a><strong><font color="#228B22">  #define MPIU_2INT       MPI_2INT</font></strong>
<a name="line1037">1037: </a><strong><font color="#228B22">  #define MPIU_INT_MPIINT MPI_2INT</font></strong>
<a name="line1038">1038: </a><font color="#A020F0">#endif</font>
<a name="line1039">1039: </a>PETSC_EXTERN MPI_Datatype MPI_4INT;
<a name="line1040">1040: </a>PETSC_EXTERN MPI_Datatype MPIU_4INT;

<a name="line1042">1042: </a><strong><font color="#4169E1"><a name="PetscPowInt"></a>static inline <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> PetscPowInt(<a href="../manualpages/Sys/PetscInt.html">PetscInt</a> base, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> power)</font></strong>
<a name="line1043">1043: </a>{
<a name="line1044">1044: </a>  <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> result = 1;
<a name="line1045">1045: </a>  <font color="#4169E1">while</font> (power) {
<a name="line1046">1046: </a>    <font color="#4169E1">if</font> (power &amp; 1) result *= base;
<a name="line1047">1047: </a>    power &gt;&gt;= 1;
<a name="line1048">1048: </a>    <font color="#4169E1">if</font> (power) base *= base;
<a name="line1049">1049: </a>  }
<a name="line1050">1050: </a>  <font color="#4169E1">return</font> result;
<a name="line1051">1051: </a>}

<a name="line1053">1053: </a><strong><font color="#4169E1"><a name="PetscPowInt64"></a>static inline PetscInt64 PetscPowInt64(<a href="../manualpages/Sys/PetscInt.html">PetscInt</a> base, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> power)</font></strong>
<a name="line1054">1054: </a>{
<a name="line1055">1055: </a>  PetscInt64 result = 1;
<a name="line1056">1056: </a>  <font color="#4169E1">while</font> (power) {
<a name="line1057">1057: </a>    <font color="#4169E1">if</font> (power &amp; 1) result *= base;
<a name="line1058">1058: </a>    power &gt;&gt;= 1;
<a name="line1059">1059: </a>    <font color="#4169E1">if</font> (power) base *= base;
<a name="line1060">1060: </a>  }
<a name="line1061">1061: </a>  <font color="#4169E1">return</font> result;
<a name="line1062">1062: </a>}

<a name="line1064">1064: </a><strong><font color="#4169E1"><a name="PetscPowRealInt"></a>static inline <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> PetscPowRealInt(<a href="../manualpages/Sys/PetscReal.html">PetscReal</a> base, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> power)</font></strong>
<a name="line1065">1065: </a>{
<a name="line1066">1066: </a>  <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> result = 1;
<a name="line1067">1067: </a>  <font color="#4169E1">if</font> (power &lt; 0) {
<a name="line1068">1068: </a>    power = -power;
<a name="line1069">1069: </a>    base  = ((<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)1) / base;
<a name="line1070">1070: </a>  }
<a name="line1071">1071: </a>  <font color="#4169E1">while</font> (power) {
<a name="line1072">1072: </a>    <font color="#4169E1">if</font> (power &amp; 1) result *= base;
<a name="line1073">1073: </a>    power &gt;&gt;= 1;
<a name="line1074">1074: </a>    <font color="#4169E1">if</font> (power) base *= base;
<a name="line1075">1075: </a>  }
<a name="line1076">1076: </a>  <font color="#4169E1">return</font> result;
<a name="line1077">1077: </a>}

<a name="line1079">1079: </a><strong><font color="#4169E1"><a name="PetscPowScalarInt"></a>static inline <a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> PetscPowScalarInt(<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> base, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> power)</font></strong>
<a name="line1080">1080: </a>{
<a name="line1081">1081: </a>  <a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> result = (<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)1;
<a name="line1082">1082: </a>  <font color="#4169E1">if</font> (power &lt; 0) {
<a name="line1083">1083: </a>    power = -power;
<a name="line1084">1084: </a>    base  = ((<a href="../manualpages/Sys/PetscReal.html">PetscReal</a>)1) / base;
<a name="line1085">1085: </a>  }
<a name="line1086">1086: </a>  <font color="#4169E1">while</font> (power) {
<a name="line1087">1087: </a>    <font color="#4169E1">if</font> (power &amp; 1) result *= base;
<a name="line1088">1088: </a>    power &gt;&gt;= 1;
<a name="line1089">1089: </a>    <font color="#4169E1">if</font> (power) base *= base;
<a name="line1090">1090: </a>  }
<a name="line1091">1091: </a>  <font color="#4169E1">return</font> result;
<a name="line1092">1092: </a>}

<a name="line1094">1094: </a><strong><font color="#4169E1"><a name="PetscPowScalarReal"></a>static inline <a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> PetscPowScalarReal(<a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> base, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> power)</font></strong>
<a name="line1095">1095: </a>{
<a name="line1096">1096: </a>  <a href="../manualpages/Sys/PetscScalar.html">PetscScalar</a> cpower = power;
<a name="line1097">1097: </a>  <font color="#4169E1">return</font> PetscPowScalar(base, cpower);
<a name="line1098">1098: </a>}

<a name="line1100">1100: </a><font color="#B22222">/*MC</font>
<a name="line1101">1101: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscApproximateLTE.html">PetscApproximateLTE</a> - Performs a less than or equal to on a given constant with a fudge for floating point numbers</font>

<a name="line1103">1103: </a><font color="#B22222">   Synopsis:</font>
<a name="line1104">1104: </a>#include <A href="../include/petscmath.h.html">&lt;petscmath.h&gt;</A>
<a name="line1105">1105: </a><font color="#B22222">   bool <a href="../manualpages/Sys/PetscApproximateLTE.html">PetscApproximateLTE</a>(<a href="../manualpages/Sys/PetscReal.html">PetscReal</a> x,constant float)</font>

<a name="line1107">1107: </a><font color="#B22222">   Not Collective</font>

<a name="line1109">1109: </a><font color="#B22222">   Input Parameters:</font>
<a name="line1110">1110: </a><font color="#B22222">+   x - the variable</font>
<a name="line1111">1111: </a><font color="#B22222">-   b - the constant float it is checking if `x` is less than or equal to</font>

<a name="line1113">1113: </a><font color="#B22222">   Level: advanced</font>

<a name="line1115">1115: </a><font color="#B22222">   Notes:</font>
<a name="line1116">1116: </a><font color="#B22222">   The fudge factor is the value `<a href="../manualpages/Sys/PETSC_SMALL.html">PETSC_SMALL</a>`</font>

<a name="line1118">1118: </a><font color="#B22222">   The constant numerical value is automatically set to the appropriate precision of PETSc so can just be provided as, for example, 3.2</font>

<a name="line1120">1120: </a><font color="#B22222">   This is used in several examples for setting initial conditions based on coordinate values that are computed with i*h that produces inexact</font>
<a name="line1121">1121: </a><font color="#B22222">   floating point results.</font>

<a name="line1123">1123: </a><font color="#B22222">   Example\:</font>
<a name="line1124">1124: </a><font color="#B22222">.vb</font>
<a name="line1125">1125: </a><font color="#B22222">  <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> x;</font>
<a name="line1126">1126: </a><font color="#B22222">  if (<a href="../manualpages/Sys/PetscApproximateLTE.html">PetscApproximateLTE</a>(x, 3.2)) { // replaces if (x &lt;= 3.2) {</font>
<a name="line1127">1127: </a><font color="#B22222">.ve</font>

<a name="line1129">1129: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscMax.html">PetscMax</a>()`, `<a href="../manualpages/Sys/PetscMin.html">PetscMin</a>()`, `<a href="../manualpages/Sys/PetscAbsInt.html">PetscAbsInt</a>()`, `<a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>()`, `<a href="../manualpages/Sys/PetscApproximateGTE.html">PetscApproximateGTE</a>()`</font>
<a name="line1130">1130: </a><font color="#B22222">M*/</font>
<a name="line1131">1131: </a><strong><font color="#228B22">#define <a href="../manualpages/Sys/PetscApproximateLTE.html">PetscApproximateLTE</a>(x, b) ((x) &lt;= (<a href="../manualpages/Sys/PetscRealConstant.html">PetscRealConstant</a>(b) + <a href="../manualpages/Sys/PETSC_SMALL.html">PETSC_SMALL</a>))</font></strong>

<a name="line1133">1133: </a><font color="#B22222">/*MC</font>
<a name="line1134">1134: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscApproximateGTE.html">PetscApproximateGTE</a> - Performs a greater than or equal to on a given constant with a fudge for floating point numbers</font>

<a name="line1136">1136: </a><font color="#B22222">   Synopsis:</font>
<a name="line1137">1137: </a>#include <A href="../include/petscmath.h.html">&lt;petscmath.h&gt;</A>
<a name="line1138">1138: </a><font color="#B22222">   bool <a href="../manualpages/Sys/PetscApproximateGTE.html">PetscApproximateGTE</a>(<a href="../manualpages/Sys/PetscReal.html">PetscReal</a> x,constant float)</font>

<a name="line1140">1140: </a><font color="#B22222">   Not Collective</font>

<a name="line1142">1142: </a><font color="#B22222">   Input Parameters:</font>
<a name="line1143">1143: </a><font color="#B22222">+   x - the variable</font>
<a name="line1144">1144: </a><font color="#B22222">-   b - the constant float it is checking if `x` is greater than or equal to</font>

<a name="line1146">1146: </a><font color="#B22222">   Level: advanced</font>

<a name="line1148">1148: </a><font color="#B22222">   Notes:</font>
<a name="line1149">1149: </a><font color="#B22222">   The fudge factor is the value `<a href="../manualpages/Sys/PETSC_SMALL.html">PETSC_SMALL</a>`</font>

<a name="line1151">1151: </a><font color="#B22222">   The constant numerical value is automatically set to the appropriate precision of PETSc so can just be provided as, for example, 3.2</font>

<a name="line1153">1153: </a><font color="#B22222">   This is used in several examples for setting initial conditions based on coordinate values that are computed with i*h that produces inexact</font>
<a name="line1154">1154: </a><font color="#B22222">   floating point results.</font>

<a name="line1156">1156: </a><font color="#B22222">   Example\:</font>
<a name="line1157">1157: </a><font color="#B22222">.vb</font>
<a name="line1158">1158: </a><font color="#B22222">  <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> x;</font>
<a name="line1159">1159: </a><font color="#B22222">  if (<a href="../manualpages/Sys/PetscApproximateGTE.html">PetscApproximateGTE</a>(x, 3.2)) {  // replaces if (x &gt;= 3.2) {</font>
<a name="line1160">1160: </a><font color="#B22222">.ve</font>

<a name="line1162">1162: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscMax.html">PetscMax</a>()`, `<a href="../manualpages/Sys/PetscMin.html">PetscMin</a>()`, `<a href="../manualpages/Sys/PetscAbsInt.html">PetscAbsInt</a>()`, `<a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>()`, `<a href="../manualpages/Sys/PetscApproximateLTE.html">PetscApproximateLTE</a>()`</font>
<a name="line1163">1163: </a><font color="#B22222">M*/</font>
<a name="line1164">1164: </a><strong><font color="#228B22">#define <a href="../manualpages/Sys/PetscApproximateGTE.html">PetscApproximateGTE</a>(x, b) ((x) &gt;= (<a href="../manualpages/Sys/PetscRealConstant.html">PetscRealConstant</a>(b) - <a href="../manualpages/Sys/PETSC_SMALL.html">PETSC_SMALL</a>))</font></strong>

<a name="line1166">1166: </a><font color="#B22222">/*@C</font>
<a name="line1167">1167: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscCeilInt.html">PetscCeilInt</a> - Returns the ceiling of the quotation of two positive integers</font>

<a name="line1169">1169: </a><font color="#B22222">   Not Collective</font>

<a name="line1171">1171: </a><font color="#B22222">   Input Parameters:</font>
<a name="line1172">1172: </a><font color="#B22222">+   x - the numerator</font>
<a name="line1173">1173: </a><font color="#B22222">-   y - the denominator</font>

<a name="line1175">1175: </a><font color="#B22222">   Level: advanced</font>

<a name="line1177">1177: </a><font color="#B22222">  Example\:</font>
<a name="line1178">1178: </a><font color="#B22222">.vb</font>
<a name="line1179">1179: </a><font color="#B22222">  <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> n = <a href="../manualpages/Sys/PetscCeilInt.html">PetscCeilInt</a>(10, 3); // n has the value of 4</font>
<a name="line1180">1180: </a><font color="#B22222">.ve</font>

<a name="line1182">1182: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscCeilInt64.html">PetscCeilInt64</a>()`, `<a href="../manualpages/Sys/PetscMax.html">PetscMax</a>()`, `<a href="../manualpages/Sys/PetscMin.html">PetscMin</a>()`, `<a href="../manualpages/Sys/PetscAbsInt.html">PetscAbsInt</a>()`, `<a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>()`, `<a href="../manualpages/Sys/PetscApproximateLTE.html">PetscApproximateLTE</a>()`</font>
<a name="line1183">1183: </a><font color="#B22222">@*/</font>
<a name="line1184">1184: </a><strong><font color="#4169E1"><a name="PetscCeilInt"></a>static inline <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> <a href="../manualpages/Sys/PetscCeilInt.html">PetscCeilInt</a>(<a href="../manualpages/Sys/PetscInt.html">PetscInt</a> x, <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> y)</font></strong>
<a name="line1185">1185: </a>{
<a name="line1186">1186: </a>  <font color="#4169E1">return</font> x / y + (x % y ? 1 : 0);
<a name="line1187">1187: </a>}

<a name="line1189">1189: </a><font color="#B22222">/*@C</font>
<a name="line1190">1190: </a><font color="#B22222">   <a href="../manualpages/Sys/PetscCeilInt64.html">PetscCeilInt64</a> - Returns the ceiling of the quotation of two positive integers</font>

<a name="line1192">1192: </a><font color="#B22222">   Not Collective</font>

<a name="line1194">1194: </a><font color="#B22222">   Input Parameters:</font>
<a name="line1195">1195: </a><font color="#B22222">+   x - the numerator</font>
<a name="line1196">1196: </a><font color="#B22222">-   y - the denominator</font>

<a name="line1198">1198: </a><font color="#B22222">   Level: advanced</font>

<a name="line1200">1200: </a><font color="#B22222">  Example\:</font>
<a name="line1201">1201: </a><font color="#B22222">.vb</font>
<a name="line1202">1202: </a><font color="#B22222">  PetscInt64 n = <a href="../manualpages/Sys/PetscCeilInt64.html">PetscCeilInt64</a>(10, 3); // n has the value of 4</font>
<a name="line1203">1203: </a><font color="#B22222">.ve</font>

<a name="line1205">1205: </a><font color="#B22222">.seealso: `<a href="../manualpages/Sys/PetscCeilInt.html">PetscCeilInt</a>()`, `<a href="../manualpages/Sys/PetscMax.html">PetscMax</a>()`, `<a href="../manualpages/Sys/PetscMin.html">PetscMin</a>()`, `<a href="../manualpages/Sys/PetscAbsInt.html">PetscAbsInt</a>()`, `<a href="../manualpages/Sys/PetscAbsReal.html">PetscAbsReal</a>()`, `<a href="../manualpages/Sys/PetscApproximateLTE.html">PetscApproximateLTE</a>()`</font>
<a name="line1206">1206: </a><font color="#B22222">@*/</font>
<a name="line1207">1207: </a><strong><font color="#4169E1"><a name="PetscCeilInt64"></a>static inline PetscInt64 <a href="../manualpages/Sys/PetscCeilInt64.html">PetscCeilInt64</a>(PetscInt64 x, PetscInt64 y)</font></strong>
<a name="line1208">1208: </a>{
<a name="line1209">1209: </a>  <font color="#4169E1">return</font> x / y + (x % y ? 1 : 0);
<a name="line1210">1210: </a>}

<a name="line1212">1212: </a><strong><font color="#4169E1">PETSC_EXTERN <a href="../manualpages/Sys/PetscErrorCode.html">PetscErrorCode</a> <a href="../manualpages/Sys/PetscLinearRegression.html">PetscLinearRegression</a>(<a href="../manualpages/Sys/PetscInt.html">PetscInt</a>, const <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>[], const <a href="../manualpages/Sys/PetscReal.html">PetscReal</a>[], <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *, <a href="../manualpages/Sys/PetscReal.html">PetscReal</a> *)</font></strong>;
</pre>
</body>

</html>