ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_mathSqrt::getSqrt Interface Reference

Generate and return the integer square root of an input non-negative integer. More...

Detailed Description

Generate and return the integer square root of an input non-negative integer.

See the documentation of pm_mathSqrt for more information on integer square root.

Parameters
[in]posint: The input scalar of type integer of kind any supported by the processor (e.g., IK, IK8, IK16, IK32, or IK64) containing the non-negative integer whose factoring is to be computed on return.
Returns
intSqrt : The output allocatable vector of the same type and kind as the input posint, containing the prime factors of the input integer posint.


Possible calling interfaces

use pm_mathSqrt, only: getSqrt
sqrt = getSqrt(posint)
Generate and return the integer square root of an input non-negative integer.
This module contains procedures and generic interfaces for computing the square root of integers.
Definition: pm_mathSqrt.F90:67
Warning
The condition 0 <= posint must hold for the corresponding input arguments.
This condition is verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
The input argument posint has the value attribute.
The computation of the integer square root using the linear search method can become extremely costly for large input posint values, typically larger than 10**5.
For more information, see the relevant benchmarks.
The pure procedure(s) documented herein become impure when the ParaMonte library is compiled with preprocessor macro CHECK_ENABLED=1.
By default, these procedures are pure in release build and impure in debug and testing builds.
Remarks
The procedures under discussion are elemental.
The procedures of this generic interface compute the integer square root method without undue overflow.
See also
getLog1p
get1mexp


Example usage

1program example
2
3 use pm_kind, only: IKG => IKS ! any integer kind is supported.
4 use pm_kind, only: SK, IK, LK, RKD
5 use pm_io, only: display_type
7 use pm_mathSqrt, only: getSqrt, linear, binary
8 use pm_err, only: setAsserted
9
10 implicit none
11
12 integer(IK) :: i
13 integer(IKG) :: posint, intSqrt
14
15 type(display_type) :: disp
16 disp = display_type(file = "main.out.F90")
17
18 do i = 1, 10
19 call disp%skip()
20 call disp%show("posint = getLogUnifRand(1, huge(1))")
21 posint = getLogUnifRand(1, huge(1))
22 call disp%show("posint")
23 call disp%show( posint )
24 call disp%show("intSqrt = getSqrt(posint, binary)")
25 intSqrt = getSqrt(posint, binary)
26 call disp%show("[intSqrt, floor(sqrt(real(posint, RKD)))]")
27 call disp%show( [intSqrt, floor(sqrt(real(posint, RKD)))] )
28 call disp%show("call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))")
29 call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
30 call disp%skip()
31 end do
32
33 do posint = huge(0_IKG), huge(0_IKG) - 10_IKG, -1_IKG
34 call disp%skip()
35 call disp%show("posint")
36 call disp%show( posint )
37 call disp%show("intSqrt = getSqrt(posint, binary)")
38 intSqrt = getSqrt(posint, binary)
39 call disp%show("[intSqrt, floor(sqrt(real(posint, RKD)))]")
40 call disp%show( [intSqrt, floor(sqrt(real(posint, RKD)))] )
41 call disp%show("call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))")
42 call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
43 call disp%skip()
44 end do
45
46 do posint = 0_IKG, 10_IKG
47 call disp%skip()
48 call disp%show("posint")
49 call disp%show( posint )
50 call disp%show("intSqrt = getSqrt(posint, binary)")
51 intSqrt = getSqrt(posint, binary)
52 call disp%show("[intSqrt, getSqrt(posint, linear), floor(sqrt(real(posint, RKD)))]")
53 call disp%show( [intSqrt, getSqrt(posint, linear), floor(sqrt(real(posint, RKD)))] )
54 call disp%show("call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))")
55 call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
56 call disp%skip()
57 end do
58
59end program example
Generate and return a scalar (or array of arbitrary rank) of random value(s) from the LogUniform dist...
Verify the input assertion holds and if it does not, print the (optional) input message on stdout and...
Definition: pm_err.F90:735
Generate and return an object of type stop_type with the user-specified input attributes.
Definition: pm_err.F90:1618
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11508
This module contains classes and procedures for computing various statistical quantities related to t...
This module contains classes and procedures for reporting and handling errors.
Definition: pm_err.F90:52
This module contains classes and procedures for input/output (IO) or generic display operations on st...
Definition: pm_io.F90:252
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
Definition: pm_io.F90:11393
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
integer, parameter IKS
The single-precision integer kind in Fortran mode. On most platforms, this is a 32-bit integer kind.
Definition: pm_kind.F90:563
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
integer, parameter RKD
The double precision real kind in Fortran mode. On most platforms, this is an 64-bit real kind.
Definition: pm_kind.F90:568
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
Generate and return an object of type display_type.
Definition: pm_io.F90:10282

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2posint = getLogUnifRand(1, huge(1))
3posint
4+96789224
5intSqrt = getSqrt(posint, binary)
6[intSqrt, floor(sqrt(real(posint, RKD)))]
7+9838, +9838
8call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
9
10
11posint = getLogUnifRand(1, huge(1))
12posint
13+858239
14intSqrt = getSqrt(posint, binary)
15[intSqrt, floor(sqrt(real(posint, RKD)))]
16+926, +926
17call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
18
19
20posint = getLogUnifRand(1, huge(1))
21posint
22+346450976
23intSqrt = getSqrt(posint, binary)
24[intSqrt, floor(sqrt(real(posint, RKD)))]
25+18613, +18613
26call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
27
28
29posint = getLogUnifRand(1, huge(1))
30posint
31+52148880
32intSqrt = getSqrt(posint, binary)
33[intSqrt, floor(sqrt(real(posint, RKD)))]
34+7221, +7221
35call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
36
37
38posint = getLogUnifRand(1, huge(1))
39posint
40+500784672
41intSqrt = getSqrt(posint, binary)
42[intSqrt, floor(sqrt(real(posint, RKD)))]
43+22378, +22378
44call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
45
46
47posint = getLogUnifRand(1, huge(1))
48posint
49+5102908
50intSqrt = getSqrt(posint, binary)
51[intSqrt, floor(sqrt(real(posint, RKD)))]
52+2258, +2258
53call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
54
55
56posint = getLogUnifRand(1, huge(1))
57posint
58+690859
59intSqrt = getSqrt(posint, binary)
60[intSqrt, floor(sqrt(real(posint, RKD)))]
61+831, +831
62call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
63
64
65posint = getLogUnifRand(1, huge(1))
66posint
67+5425
68intSqrt = getSqrt(posint, binary)
69[intSqrt, floor(sqrt(real(posint, RKD)))]
70+73, +73
71call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
72
73
74posint = getLogUnifRand(1, huge(1))
75posint
76+35
77intSqrt = getSqrt(posint, binary)
78[intSqrt, floor(sqrt(real(posint, RKD)))]
79+5, +5
80call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
81
82
83posint = getLogUnifRand(1, huge(1))
84posint
85+112842752
86intSqrt = getSqrt(posint, binary)
87[intSqrt, floor(sqrt(real(posint, RKD)))]
88+10622, +10622
89call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
90
91
92posint
93+2147483647
94intSqrt = getSqrt(posint, binary)
95[intSqrt, floor(sqrt(real(posint, RKD)))]
96+46340, +46340
97call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
98
99
100posint
101+2147483646
102intSqrt = getSqrt(posint, binary)
103[intSqrt, floor(sqrt(real(posint, RKD)))]
104+46340, +46340
105call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
106
107
108posint
109+2147483645
110intSqrt = getSqrt(posint, binary)
111[intSqrt, floor(sqrt(real(posint, RKD)))]
112+46340, +46340
113call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
114
115
116posint
117+2147483644
118intSqrt = getSqrt(posint, binary)
119[intSqrt, floor(sqrt(real(posint, RKD)))]
120+46340, +46340
121call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
122
123
124posint
125+2147483643
126intSqrt = getSqrt(posint, binary)
127[intSqrt, floor(sqrt(real(posint, RKD)))]
128+46340, +46340
129call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
130
131
132posint
133+2147483642
134intSqrt = getSqrt(posint, binary)
135[intSqrt, floor(sqrt(real(posint, RKD)))]
136+46340, +46340
137call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
138
139
140posint
141+2147483641
142intSqrt = getSqrt(posint, binary)
143[intSqrt, floor(sqrt(real(posint, RKD)))]
144+46340, +46340
145call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
146
147
148posint
149+2147483640
150intSqrt = getSqrt(posint, binary)
151[intSqrt, floor(sqrt(real(posint, RKD)))]
152+46340, +46340
153call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
154
155
156posint
157+2147483639
158intSqrt = getSqrt(posint, binary)
159[intSqrt, floor(sqrt(real(posint, RKD)))]
160+46340, +46340
161call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
162
163
164posint
165+2147483638
166intSqrt = getSqrt(posint, binary)
167[intSqrt, floor(sqrt(real(posint, RKD)))]
168+46340, +46340
169call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
170
171
172posint
173+2147483637
174intSqrt = getSqrt(posint, binary)
175[intSqrt, floor(sqrt(real(posint, RKD)))]
176+46340, +46340
177call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
178
179
180posint
181+0
182intSqrt = getSqrt(posint, binary)
183[intSqrt, getSqrt(posint, linear), floor(sqrt(real(posint, RKD)))]
184+0, +0, +0
185call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
186
187
188posint
189+1
190intSqrt = getSqrt(posint, binary)
191[intSqrt, getSqrt(posint, linear), floor(sqrt(real(posint, RKD)))]
192+1, +1, +1
193call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
194
195
196posint
197+2
198intSqrt = getSqrt(posint, binary)
199[intSqrt, getSqrt(posint, linear), floor(sqrt(real(posint, RKD)))]
200+1, +1, +1
201call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
202
203
204posint
205+3
206intSqrt = getSqrt(posint, binary)
207[intSqrt, getSqrt(posint, linear), floor(sqrt(real(posint, RKD)))]
208+1, +1, +1
209call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
210
211
212posint
213+4
214intSqrt = getSqrt(posint, binary)
215[intSqrt, getSqrt(posint, linear), floor(sqrt(real(posint, RKD)))]
216+2, +2, +2
217call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
218
219
220posint
221+5
222intSqrt = getSqrt(posint, binary)
223[intSqrt, getSqrt(posint, linear), floor(sqrt(real(posint, RKD)))]
224+2, +2, +2
225call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
226
227
228posint
229+6
230intSqrt = getSqrt(posint, binary)
231[intSqrt, getSqrt(posint, linear), floor(sqrt(real(posint, RKD)))]
232+2, +2, +2
233call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
234
235
236posint
237+7
238intSqrt = getSqrt(posint, binary)
239[intSqrt, getSqrt(posint, linear), floor(sqrt(real(posint, RKD)))]
240+2, +2, +2
241call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
242
243
244posint
245+8
246intSqrt = getSqrt(posint, binary)
247[intSqrt, getSqrt(posint, linear), floor(sqrt(real(posint, RKD)))]
248+2, +2, +2
249call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
250
251
252posint
253+9
254intSqrt = getSqrt(posint, binary)
255[intSqrt, getSqrt(posint, linear), floor(sqrt(real(posint, RKD)))]
256+3, +3, +3
257call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
258
259
260posint
261+10
262intSqrt = getSqrt(posint, binary)
263[intSqrt, getSqrt(posint, linear), floor(sqrt(real(posint, RKD)))]
264+3, +3, +3
265call setAsserted(intSqrt == floor(sqrt(real(posint, RKD))))
266
267
Test:
test_pm_mathIntSqrt
Todo:
Low Priority: A binary-representation calculation of the integer square root should be added in the future.


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, April 23, 2017, 1:36 AM, Institute for Computational Engineering and Sciences (ICES), University of Texas at Austin

Definition at line 141 of file pm_mathSqrt.F90.


The documentation for this interface was generated from the following file: