Python* API Reference for Intel® Data Analytics Acceleration Library 2018 Update 1

qr_dense_distr.py

1 #===============================================================================
2 # Copyright 2014-2017 Intel Corporation
3 # All Rights Reserved.
4 #
5 # If this software was obtained under the Intel Simplified Software License,
6 # the following terms apply:
7 #
8 # The source code, information and material ("Material") contained herein is
9 # owned by Intel Corporation or its suppliers or licensors, and title to such
10 # Material remains with Intel Corporation or its suppliers or licensors. The
11 # Material contains proprietary information of Intel or its suppliers and
12 # licensors. The Material is protected by worldwide copyright laws and treaty
13 # provisions. No part of the Material may be used, copied, reproduced,
14 # modified, published, uploaded, posted, transmitted, distributed or disclosed
15 # in any way without Intel's prior express written permission. No license under
16 # any patent, copyright or other intellectual property rights in the Material
17 # is granted to or conferred upon you, either expressly, by implication,
18 # inducement, estoppel or otherwise. Any license under such intellectual
19 # property rights must be express and approved by Intel in writing.
20 #
21 # Unless otherwise agreed by Intel in writing, you may not remove or alter this
22 # notice or any other notice embedded in Materials by Intel or Intel's
23 # suppliers or licensors in any way.
24 #
25 #
26 # If this software was obtained under the Apache License, Version 2.0 (the
27 # "License"), the following terms apply:
28 #
29 # You may not use this file except in compliance with the License. You may
30 # obtain a copy of the License at http://www.apache.org/licenses/LICENSE-2.0
31 #
32 #
33 # Unless required by applicable law or agreed to in writing, software
34 # distributed under the License is distributed on an "AS IS" BASIS, WITHOUT
35 # WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
36 #
37 # See the License for the specific language governing permissions and
38 # limitations under the License.
39 #===============================================================================
40 
41 
42 
43 
44 import os
45 import sys
46 
47 from daal import step1Local, step2Master, step3Local
48 from daal.algorithms import qr
49 from daal.data_management import FileDataSource, DataSourceIface
50 
51 utils_folder = os.path.realpath(os.path.abspath(os.path.dirname(os.path.dirname(__file__))))
52 if utils_folder not in sys.path:
53  sys.path.insert(0, utils_folder)
54 from utils import printNumericTable
55 
56 DAAL_PREFIX = os.path.join('..', 'data')
57 
58 # Input data set parameters
59 nBlocks = 4
60 
61 datasetFileNames = [
62  os.path.join(DAAL_PREFIX, 'distributed', 'qr_1.csv'),
63  os.path.join(DAAL_PREFIX, 'distributed', 'qr_2.csv'),
64  os.path.join(DAAL_PREFIX, 'distributed', 'qr_3.csv'),
65  os.path.join(DAAL_PREFIX, 'distributed', 'qr_4.csv')
66 ]
67 
68 dataFromStep1ForStep2 = [0] * nBlocks
69 dataFromStep1ForStep3 = [0] * nBlocks
70 dataFromStep2ForStep3 = [0] * nBlocks
71 R = None
72 Qi = [0] * nBlocks
73 
74 
75 def computestep1Local(block):
76  global dataFromStep1ForStep2, dataFromStep1ForStep3
77 
78  # Initialize FileDataSource<CSVFeatureManager> to retrieve the input data from a .csv file
79  dataSource = FileDataSource(
80  datasetFileNames[block],
81  DataSourceIface.doAllocateNumericTable,
82  DataSourceIface.doDictionaryFromContext
83  )
84 
85  # Retrieve the input data
86  dataSource.loadDataBlock()
87 
88  # Create an algorithm to compute QR decomposition on the local node
89  algorithm = qr.Distributed(step1Local)
90 
91  algorithm.input.set(qr.data, dataSource.getNumericTable())
92 
93  # Compute QR decomposition and get OnlinePartialResult class from daal.algorithms.qr
94  pres = algorithm.compute()
95 
96  dataFromStep1ForStep2[block] = pres.get(qr.outputOfStep1ForStep2)
97  dataFromStep1ForStep3[block] = pres.get(qr.outputOfStep1ForStep3)
98 
99 
100 def computeOnMasterNode():
101  global R, dataFromStep2ForStep3
102 
103  # Create an algorithm to compute QR decomposition on the master node
104  algorithm = qr.Distributed(step2Master)
105 
106  for i in range(nBlocks):
107  algorithm.input.add(qr.inputOfStep2FromStep1, i, dataFromStep1ForStep2[i])
108 
109  # Compute QR decomposition and get DistributedPartialResult class from daal.algorithms.qr
110  pres = algorithm.compute()
111 
112  for i in range(nBlocks):
113  dataFromStep2ForStep3[i] = pres.getCollection(qr.outputOfStep2ForStep3, i)
114 
115  res = algorithm.finalizeCompute()
116  R = res.get(qr.matrixR)
117 
118 
119 def finalizeComputestep1Local(block):
120  global Qi
121 
122  # Create an algorithm to compute QR decomposition on the master node
123  algorithm = qr.Distributed(step3Local)
124 
125  algorithm.input.set(qr.inputOfStep3FromStep1, dataFromStep1ForStep3[block])
126  algorithm.input.set(qr.inputOfStep3FromStep2, dataFromStep2ForStep3[block])
127 
128  # Compute QR decomposition
129  algorithm.compute()
130 
131  res = algorithm.finalizeCompute()
132 
133  Qi[block] = res.get(qr.matrixQ)
134 
135 if __name__ == "__main__":
136 
137  for i in range(nBlocks):
138  computestep1Local(i)
139 
140  computeOnMasterNode()
141 
142  for i in range(nBlocks):
143  finalizeComputestep1Local(i)
144 
145  # Print the results
146  printNumericTable(Qi[0], "Part of orthogonal matrix Q from 1st node:", 10)
147  printNumericTable(R, "Triangular matrix R:")

For more complete information about compiler optimizations, see our Optimization Notice.