File: testRefinement.py

package info (click to toggle)
tasmanian 8.2-1
  • links: PTS, VCS
  • area: main
  • in suites: sid
  • size: 4,852 kB
  • sloc: cpp: 34,523; python: 7,039; f90: 5,080; makefile: 224; sh: 64; ansic: 8
file content (259 lines) | stat: -rw-r--r-- 12,423 bytes parent folder | download | duplicates (3)
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
import unittest
import TasmanianSG
import numpy as np

from random import shuffle

import testCommon

ttc = testCommon.TestTasCommon()

class TestTasClass(unittest.TestCase):
    '''
    Test the refinement capabilities:
    * set different refinements
    * estimate anisotropic coefficients
    * read/write refinements
    '''
    def __init__(self):
        unittest.TestCase.__init__(self, "testNothing")

    def testNothing(self):
        pass

    def checkSetClear(self):
        '''
        Set refinement and clear refinement
        '''
        grid = TasmanianSG.TasmanianSparseGrid()
        grid.makeLocalPolynomialGrid(2, 1, 2, 1, 'semi-localp')
        ttc.loadExpN2(grid)
        self.assertEqual(grid.getNumNeeded(), 0, "num needed")
        grid.setSurplusRefinement(0.0001, 0, 'classic')
        self.assertTrue((grid.getNumNeeded() > 0), "num needed")
        grid.clearRefinement()
        self.assertEqual(grid.getNumNeeded(), 0, "num needed")

    def checkAnisoCoeff(self):
        '''
        Check anisotropic coefficients
        '''
        grid = TasmanianSG.TasmanianSparseGrid()
        grid.makeGlobalGrid(2, 1, 9, 'level', 'rleja')
        aP = grid.getPoints()
        aV = np.exp(aP[:,0] + aP[:,1]**2)
        grid.loadNeededPoints(aV.reshape([aP.shape[0], 1]))
        aC = grid.estimateAnisotropicCoefficients('iptotal', 0)
        self.assertEqual(len(aC.shape), 1, 'dimensions of the estimated anisotropic weight')
        self.assertEqual(aC.shape[0], 2, 'dimensions of the estimated anisotropic weight')
        self.assertLess(np.abs(float(aC[0]) / float(aC[1]) - 2.0), 0.2, 'wrong anisotropic weights estimated')
        aC = grid.estimateAnisotropicCoefficients('ipcurved', 0)
        self.assertEqual(len(aC.shape), 1, 'dimensions of the estimated anisotropic weight')
        self.assertEqual(aC.shape[0], 4, 'dimensions of the estimated anisotropic weight')
        self.assertLess(np.abs(float(aC[0]) / float(aC[1]) - 2.0), 0.2, 'wrong anisotropic weights estimated, alpha curved')
        self.assertLess(aC[2], 0.0, 'wrong anisotropic weights estimated, beta 1')
        self.assertLess(aC[3], 0.0, 'wrong anisotropic weights estimated, beta 2')

    def checkLocalpSurplus(self):
        '''
        Check surplus refinement for local polynomial grids
        '''
        grid = TasmanianSG.TasmanianSparseGrid()
        grid.makeLocalPolynomialGrid(2, 1, 4, 1, 'semi-localp')
        ttc.loadExpN2(grid)
        aPoints = grid.getPoints()
        aScale = np.array([[1.0 if aPoints[i,0] > 0.0 else 0.0] for i in range(aPoints.shape[0])])
        grid.setSurplusRefinement(1.E-9, 0, 'classic', [], aScale)
        aNeeded = grid.getNeededPoints()
        for iI in range(aNeeded.shape[0]):
            self.assertLess(0.0, aNeeded[iI, 0], 'wrong set of needed points after rescaling')


    def checkFileIO(self):
        '''
        Read/Write regular refinement.
        '''
        grid = TasmanianSG.TasmanianSparseGrid()

        sRefinementIOGrids = ["grid.makeGlobalGrid(2, 1, 2, 'level', 'clenshaw-curtis')",
                              "grid.makeSequenceGrid(2, 1, 2, 'level', 'rleja')"]
        for sTest in sRefinementIOGrids:
            exec(sTest)
            ttc.loadExpN2(grid)
            grid.setAnisotropicRefinement('iptotal', 20, 0)
            self.assertTrue(grid.getNumNeeded() > 0, 'did not refine')
            gridB = TasmanianSG.TasmanianSparseGrid()
            gridB.copyGrid(grid)
            ttc.compareGrids(grid, gridB)

            grid.write("refTestFlename.grid", bUseBinaryFormat = True)
            gridB.makeLocalPolynomialGrid(1, 1, 0, 1)
            gridB.read("refTestFlename.grid")
            ttc.compareGrids(grid, gridB)

            grid.write("refTestFlename.grid", bUseBinaryFormat = False)
            gridB.makeLocalPolynomialGrid(1, 1, 0, 1)
            gridB.read("refTestFlename.grid")
            ttc.compareGrids(grid, gridB)

    def checkConstruction(self):
        '''
        Test read/write when using construction.
        '''
        llTest = ["gridA.makeGlobalGrid(3, 2, 2, 'level', 'clenshaw-curtis'); gridB.makeGlobalGrid(3, 2, 2, 'level', 'clenshaw-curtis')",
                  "gridA.makeSequenceGrid(3, 2, 4, 'level', 'leja'); gridB.makeSequenceGrid(3, 2, 4, 'level', 'leja')",
                  "gridA.makeLocalPolynomialGrid(3, 2, 2); gridB.makeLocalPolynomialGrid(3, 2, 2)",
                  "gridA.makeWaveletGrid(3, 2, 2); gridB.makeWaveletGrid(3, 2, 2)",
                  "gridA.makeFourierGrid(3, 2, 2, 'level'); gridB.makeFourierGrid(3, 2, 2, 'level')",]

        for sMakeGrids in llTest:
            for sFormat in [False, True]: # test binary and ascii format
                gridA = TasmanianSG.TasmanianSparseGrid()
                gridB = TasmanianSG.TasmanianSparseGrid()
                gridC = TasmanianSG.TasmanianSparseGrid()

                exec(sMakeGrids)

                gridA.beginConstruction()
                gridB.beginConstruction()
                #gridA.printStats()

                gridB.write("testSave", bUseBinaryFormat = sFormat)
                gridB.makeSequenceGrid(1, 1, 0, "level", "rleja") # clean the grid
                gridB.read("testSave")
                ttc.compareGrids(gridA, gridB)
                gridC.copyGrid(gridA)
                ttc.compareGrids(gridA, gridC)

                for t in range(5): # use 5 iterations
                    if (gridA.isLocalPolynomial() or gridA.isWavelet()):
                        aPointsA = gridA.getCandidateConstructionPointsSurplus(1.E-4, "fds")
                        aPointsB = gridB.getCandidateConstructionPointsSurplus(1.E-4, "fds")
                        aPointsC = gridC.getCandidateConstructionPointsSurplus(1.E-4, "fds")
                    else:
                        aPointsA = gridA.getCandidateConstructionPoints("level", 0)
                        aPointsB = gridB.getCandidateConstructionPoints("level", 0)
                        aPointsC = gridC.getCandidateConstructionPoints("level", 0)
                    np.testing.assert_almost_equal(aPointsA, aPointsB, decimal=11)
                    np.testing.assert_almost_equal(aPointsA, aPointsC, decimal=11)

                    iNumPoints = int(aPointsA.shape[0] / 2)
                    if (iNumPoints > 32): iNumPoints = 32

                    # use the first samples (up to 32) and shuffle the order
                    # add one of the samples further in the list
                    liSamples = list(range(iNumPoints + 1))
                    shuffle(liSamples)
                    for iI in range(len(liSamples)):
                        if (liSamples[iI] == iNumPoints):
                            liSamples[iI] = iNumPoints + 1
                    #liSamples = map(lambda i: i if i < iNumPoints else iNumPoints + 1, liSamples)

                    for iI in liSamples: # compute and load the samples
                        aPoint = aPointsA[iI, :]
                        aValue = np.array([np.exp(aPoint[0] + aPoint[1]), 1.0 / ((aPoint[0] - 1.3) * (aPoint[1] - 1.6) * (aPoint[2] - 2.0))])

                        gridA.loadConstructedPoint(aPoint, aValue)
                        gridB.loadConstructedPoint(aPoint, aValue)
                        gridC.loadConstructedPoint(aPoint, aValue)

                    # using straight construction or read/write should produce the same result
                    ttc.compareGrids(gridA, gridC)
                    gridB.write("testSave", bUseBinaryFormat = sFormat)
                    gridB.makeSequenceGrid(1, 1, 0, "level", "rleja")
                    gridB.read("testSave")
                    ttc.compareGrids(gridA, gridB)
                    gridC.copyGrid(gridA)
                    ttc.compareGrids(gridA, gridC)

                gridA.finishConstruction()
                gridB.finishConstruction()

                gridB.write("testSave", bUseBinaryFormat = sFormat)
                gridB.makeSequenceGrid(1, 1, 0, "level", "rleja")
                gridB.read("testSave")
                ttc.compareGrids(gridA, gridB)
                gridC.copyGrid(gridA)
                ttc.compareGrids(gridA, gridC)

        # check multi-point load
        gridA = TasmanianSG.TasmanianSparseGrid()
        gridA.makeLocalPolynomialGrid(3, 2, 4);
        ttc.loadExpN2(gridA)

        gridB = TasmanianSG.TasmanianSparseGrid()
        gridB.makeLocalPolynomialGrid(3, 2, 0)

        gridB.beginConstruction()
        aX = gridA.getPoints()
        aY = gridA.evaluateBatch(aX)
        gridB.loadConstructedPoint(aX, aY)
        gridB.finishConstruction()
        ttc.compareGrids(gridA, gridB)

        # check some mem-leaks and crashes (correctness is elsewhere)
        gridA = TasmanianSG.TasmanianSparseGrid()
        gridA.makeLocalPolynomialGrid(2, 5, 0)
        gridA.beginConstruction()
        gridA.loadConstructedPoint(np.empty([0, 2]), np.empty([0, 5])) # empty input, check for crash

        gridA.makeLocalPolynomialGrid(2, 1, 1)
        gridA.loadNeededPoints(np.ones([5, 1]))
        gridA.beginConstruction()
        aPoints = gridA.getCandidateConstructionPointsSurplus(1.E-4, "classic") # should generate empty output
        np.testing.assert_almost_equal(aPoints, np.empty([0, 0]), 14, "failed to generate empty list of construction points", True)

        gridA.makeLocalPolynomialGrid(2, 1, 0)
        gridA.loadNeededPoints(np.ones([1, 1]))
        gridA.beginConstruction()
        aPoints = gridA.getCandidateConstructionPointsSurplus(1.E-4, "classic", 0, [], np.array([[1.E-6]])) # should generate empty output
        np.testing.assert_almost_equal(aPoints, np.empty([0, 0]), 14, "failed to generate empty list of construction points", True)

        gridA.makeGlobalGrid(2, 1, 1, "tensor", "clenshaw-curtis")
        gridA.loadNeededPoints(np.ones([9, 1]))
        gridA.beginConstruction()
        aPoints = gridA.getCandidateConstructionPoints("ipcurved", [5, 5, 2, 2], [1, 1]) # should generate empty output
        np.testing.assert_almost_equal(aPoints, np.empty([0, 0]), 14, "failed to generate empty list of construction points", True)

    def checkRemovePoints(self):
        '''
        tests removePointsByHierarchicalCoefficient()
        '''
        grid = TasmanianSG.makeLocalPolynomialGrid(2, 1, 1)
        aPoints = grid.getNeededPoints()
        grid.loadNeededValues(np.exp(-aPoints[:,0]**2 -0.5*aPoints[:,1]**2).reshape((grid.getNumNeeded(), 1)))

        reduced = TasmanianSG.TasmanianSparseGrid()

        reduced.copyGrid(grid)
        reduced.removePointsByHierarchicalCoefficient(0.6)
        self.assertEqual(reduced.getNumPoints(), 3, "failed to remove points with threshold 0.6")
        np.testing.assert_almost_equal(reduced.getLoadedPoints(), np.array([[0.0, 0.0], [-1.0, 0.0], [1.0, 0.0]]), 14, "failed reduce 1", True)

        reduced.copyGrid(grid)
        reduced.removePointsByHierarchicalCoefficient(0.7)
        self.assertEqual(reduced.getNumPoints(), 1, "failed to remove points with threshold 0.7")
        np.testing.assert_almost_equal(reduced.getLoadedPoints(), np.array([[0.0, 0.0]]), 14, "failed reduce 2", True)

        reduced.copyGrid(grid)
        reduced.removePointsByHierarchicalCoefficient(0.0, iNumKeep = 3)
        self.assertEqual(reduced.getNumPoints(), 3, "failed to remove points down to 3")
        np.testing.assert_almost_equal(reduced.getLoadedPoints(), np.array([[0.0, 0.0], [-1.0, 0.0], [1.0, 0.0]]), 14, "failed reduce 3", True)

        reduced.copyGrid(grid)
        reduced.removePointsByHierarchicalCoefficient(0.0, iNumKeep = 1)
        self.assertEqual(reduced.getNumPoints(), 1, "failed to remove points down to 1")
        np.testing.assert_almost_equal(reduced.getLoadedPoints(), np.array([[0.0, 0.0]]), 14, "failed reduce 4", True)

        reduced.copyGrid(grid)
        reduced.removePointsByHierarchicalCoefficient(0.0, aScaleCorrection = np.array([[1.0], [1.0], [1.0], [0.1], [0.1]]), iNumKeep = 3)
        self.assertEqual(reduced.getNumPoints(), 3, "failed to remove corrected points")
        np.testing.assert_almost_equal(reduced.getLoadedPoints(), np.array([[0.0, 0.0], [0.0, -1.0], [0.0, 1.0]]), 14, "failed reduce 5", True)

    def performRefinementTest(self):
        self.checkSetClear()
        self.checkAnisoCoeff()
        self.checkLocalpSurplus()
        self.checkFileIO()
        self.checkConstruction()
        self.checkRemovePoints()