/
FunctionalRayleighQuotientSeparable.py
390 lines (262 loc) · 13.9 KB
/
FunctionalRayleighQuotientSeparable.py
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
import FEM
import math
import scipy as np
import sys
#
# Fem Functional : most imp function in this code
#
class FEMFunctional:
def __init__(self,
fem,
rankVeff,
sigma,
tuckerDecomposedVeff):
self.fem = fem
self.rankVeff = rankVeff
self.sigma = sigma
self.umat = tuckerDecomposedVeff[0]
self.vmat = tuckerDecomposedVeff[1]
self.wmat = tuckerDecomposedVeff[2]
self.vShapeFunctionAtQuadPoints = np.tile(fem.shapeFunctionAtQuadPoints,(1,fem.numberElements))
self.vShapeFunctionDerivsAtQuadPoints = np.tile(fem.shapeFunctionDerivsAtQuadPoints,(1,fem.numberElements))
# extra data structs used while computing force
numberNodes = fem.numberNodes
numberNodesPerElement = fem.getNumberNodesPerElement()
numberQuadraturePoints = fem.getNumberQuadPointsPerElement()
numberElements = fem.numberElements
self.Fx = np.zeros((numberNodes))
self.Fy = np.zeros((numberNodes))
self.Fz = np.zeros((numberNodes))
self.X0 = np.zeros((numberNodesPerElement,numberQuadraturePoints*numberElements))
self.Y0 = np.zeros((numberNodesPerElement,numberQuadraturePoints*numberElements))
self.Z0 = np.zeros((numberNodesPerElement,numberQuadraturePoints*numberElements))
#
# static condensation
#
def condenseMatrix(self,H):
# applyBoundaryConditions on Hx Hy Hz
H = np.delete(H,0,0)
H = np.delete(H,-1,0)
H = np.delete(H,0,1)
H = np.delete(H,-1,1)
return H
#
#
#
def computeIntegralPsiSquare(self,
psiQuadValues,
DPsiQuadValues):
fem = self.fem
#
#evaluate normPsix
#
quadValuesx = psiQuadValues[0,:]
normPsix = fem.getIntegral1D(quadValuesx*quadValuesx)
#
#evaluate normPsiy
#
quadValuesy = psiQuadValues[1,:]
normPsiy = fem.getIntegral1D(quadValuesy*quadValuesy)
#
#evaluate normPsiz
#
quadValuesz= psiQuadValues[2,:]
normPsiz = fem.getIntegral1D(quadValuesz*quadValuesz)
#
#evaluate normDPsix
#
quadValuesDpsix = DPsiQuadValues[0,:]
normDPsix = fem.getInvIntegral1D(quadValuesDpsix*quadValuesDpsix)
#
#evaluate normDPsiy
#
quadValuesDpsiy = DPsiQuadValues[1,:]
normDPsiy = fem.getInvIntegral1D(quadValuesDpsiy*quadValuesDpsiy)
#
#evaluate normDPsiz
#
quadValuesDpsiz = DPsiQuadValues[2,:]
normDPsiz = fem.getInvIntegral1D(quadValuesDpsiz*quadValuesDpsiz)
return (normPsix,normPsiy,normPsiz,normDPsix,normDPsiy,normDPsiz)
def computeSeparablePotentialsUsingTuckerVeff(self,
psiQuadValues):
fem = self.fem
numberQuadPoints = fem.getNumberQuadPointsPerElement()
numberElements = fem.numberElements
rankVeff = self.rankVeff
umat = self.umat
vmat = self.vmat
wmat = self.wmat
sigma = self.sigma
#
#allocate storage for some integrals
#
intUiPsixPsix = np.zeros((rankVeff))
intVjPsiyPsiy = np.zeros((rankVeff))
intWkPsizPsiz = np.zeros((rankVeff))
#
#first precompute integrals of u_i(x)*psix*psix where i = 1 to R
#R denotes rank used in tucker decomposition of Veff
#
a = psiQuadValues[0,:]
b = psiQuadValues[1,:]
c = psiQuadValues[2,:]
aa = a*a
bb = b*b
cc = c*c
for i in range(0,rankVeff):
quadValuesx = aa*umat[:,i]
intUiPsixPsix[i] = fem.getIntegral1D(quadValuesx)
#
#second precompute integrals of v_j(x)*psiy*psiy where j = 1 to R
#R denotes rank used in tucker decomposition of Veff
#
for j in range(0,rankVeff):
quadValuesy = bb*vmat[:,j]
intVjPsiyPsiy[j] = fem.getIntegral1D(quadValuesy)
#
#second precompute integrals of w_k(x)*psiz*psiz where k = 1 to R
#R denotes rank used in tucker decomposition of Veff
#
for k in range(0,rankVeff):
quadValuesz = cc*wmat[:,k]
intWkPsizPsiz[k] = fem.getIntegral1D(quadValuesz)
vx = np.zeros((numberElements*numberQuadPoints))
vy = np.zeros((numberElements*numberQuadPoints))
vz = np.zeros((numberElements*numberQuadPoints))
# e.g for a quadratic potential v = x^2+y^2+z^2
#vx = 0.5*x*x*integral(psiy*psiy)*integral(psiz*psiz) +
# 0.5*integral(y*y*psiy*psiy)*integral(psiz*psiz) +
# 0.5*integral(psiy*psiy)*integral(z*z*psiz*psiz)
#vy and vz similarly follows
for index in range(0,numberElements*numberQuadPoints):
dyadX = np.outer(np.outer(umat[index,:],intVjPsiyPsiy),intWkPsizPsiz)
vx[index]+= np.sum(sigma*dyadX)
dyadY = np.outer(np.outer(vmat[index,:],intUiPsixPsix),intWkPsizPsiz)
vy[index] += np.sum(sigma*dyadY)
dyadZ = np.outer(np.outer(wmat[index,:],intUiPsixPsix),intVjPsiyPsiy)
vz[index] += np.sum(sigma*dyadZ)
return(vx,vy,vz)
def computeVectorizedForce(self,lagrangeMultiplier,nodalFields):
(psiQuadValues,DPsiQuadValues) = self.fem.computeFieldsAtAllQuadPoints(nodalFields)
norms = self.computeIntegralPsiSquare(psiQuadValues,DPsiQuadValues)
(vx,vy,vz) = self.computeSeparablePotentialsUsingTuckerVeff(psiQuadValues)
# data members
fem = self.fem
numberNodes = fem.getNumberNodes()
numberNodesPerElement = fem.getNumberNodesPerElement()
numberQuadraturePoints = fem.getNumberQuadPointsPerElement()
numberElements = fem.numberElements
weightQuadPointValues = fem.weightQuadPointValues
jacobianQuadPointValues = fem.jacobianQuadPointValues
invJacobianQuadPointValues = fem.invJacobianQuadPointValues
vShapeFunctionAtQuadPoints = self.vShapeFunctionAtQuadPoints
vShapeFunctionDerivsAtQuadPoints = self.vShapeFunctionDerivsAtQuadPoints
elementConnectivity = fem.elementConnectivity
Fx = self.Fx; Fy = self.Fy; Fz = self.Fz
X0 = self.X0; Y0 = self.Y0; Z0 = self.Z0
# Fx Fy Fz = 1*numNodes
Fx.fill(0)
Fy.fill(0)
Fz.fill(0)
X0.fill(0)
Y0.fill(0)
Z0.fill(0)
normPsix = norms[0]; normPsiy = norms[1]; normPsiz = norms[2];
normDPsix = norms[3]; normDPsiy = norms[4]; normDPsiz = norms[5];
cx = 1.0/(normPsiy*normPsiz)
cy = 1.0/(normPsiz*normPsix)
cz = 1.0/(normPsix*normPsiy)
Tx = 0.5*(normDPsiy/normPsiy + normDPsiz/normPsiz)
Ty = 0.5*(normDPsiz/normPsiz + normDPsix/normPsix)
Tz = 0.5*(normDPsix/normPsix + normDPsiy/normPsiy)
for iNode in xrange(0,numberNodesPerElement):
kinetic = 0.5*DPsiQuadValues[0,:]
other = cx*psiQuadValues[0,:]*vx+ (lagrangeMultiplier+Tx)*psiQuadValues[0,:]
X0[iNode,:] = weightQuadPointValues*(kinetic*invJacobianQuadPointValues*vShapeFunctionDerivsAtQuadPoints[iNode,:] + jacobianQuadPointValues*vShapeFunctionAtQuadPoints[iNode,:]*other)
kinetic = 0.5*DPsiQuadValues[1,:]
other = cy*psiQuadValues[1,:]*vy+ (lagrangeMultiplier+Ty)*psiQuadValues[1,:]
Y0[iNode,:] = weightQuadPointValues*(kinetic*invJacobianQuadPointValues*vShapeFunctionDerivsAtQuadPoints[iNode,:] + jacobianQuadPointValues*vShapeFunctionAtQuadPoints[iNode,:]*other)
kinetic = 0.5*DPsiQuadValues[2,:]
other = cz*psiQuadValues[2,:]*vz+ (lagrangeMultiplier+Tz)*psiQuadValues[2,:]
Z0[iNode,:] = weightQuadPointValues*(kinetic*invJacobianQuadPointValues*vShapeFunctionDerivsAtQuadPoints[iNode,:] + jacobianQuadPointValues*vShapeFunctionAtQuadPoints[iNode,:]*other)
for e in range(0,numberElements):
start = e*numberQuadraturePoints
end = (e+1)*numberQuadraturePoints
for iNode in xrange(0,numberNodesPerElement):
iGlobal = elementConnectivity[iNode,e]
Fx[iGlobal] += np.sum(X0[iNode,start:end])
Fy[iGlobal] += np.sum(Y0[iNode,start:end])
Fz[iGlobal] += np.sum(Z0[iNode,start:end])
F = np.concatenate((Fx[1:-1],Fy[1:-1],Fz[1:-1]))
# append constraint force
constraintForce = normPsix*normPsiy*normPsiz - 1.0
F = np.append(F,constraintForce)
return F
def generateHamiltonianGenericPotential(self,
nodalFields):
(psiQuadValues,DPsiQuadValues) = self.fem.computeFieldsAtAllQuadPoints(nodalFields)
norms = self.computeIntegralPsiSquare(psiQuadValues,DPsiQuadValues)
(vx,vy,vz) = self.computeSeparablePotentialsUsingTuckerVeff(psiQuadValues)
normx = norms[0]; normy = norms[1]; normz = norms[2];
# data members
fem = self.fem
numberNodes = fem.getNumberNodes()
numberNodesPerElement = fem.getNumberNodesPerElement()
numberQuadraturePoints = fem.getNumberQuadPointsPerElement()
numberElements = fem.numberElements
weightQuadPointValues = fem.weightQuadPointValues
jacobianQuadPointValues = fem.jacobianQuadPointValues
invJacobianQuadPointValues = fem.invJacobianQuadPointValues
shapeFunctionAtQuadPoints = fem.shapeFunctionAtQuadPoints
shapeFunctionDerivsAtQuadPoints = fem.shapeFunctionDerivsAtQuadPoints
elementConnectivity = fem.elementConnectivity
Hx = np.zeros((numberNodes,numberNodes));
Hy = Hx.copy();
Hz = Hy.copy();
M = Hx.copy();
elementStiffnessMatX = np.zeros((numberNodesPerElement,numberNodesPerElement));
elementStiffnessMatY = elementStiffnessMatX.copy()
elementStiffnessMatZ = elementStiffnessMatX.copy()
elementMassMat = elementStiffnessMatX.copy()
#
for e in range(numberElements):
start = e*numberQuadraturePoints
end = (e+1)*numberQuadraturePoints
for iNode in range(numberNodesPerElement):
shapeFunctionGradientQuadPointValuesI = shapeFunctionDerivsAtQuadPoints[iNode,:]
shapeFunctionQuadPointValuesI = shapeFunctionAtQuadPoints[iNode,:]
for jNode in range(numberNodesPerElement):
shapeFunctionGradientQuadPointValuesJ = shapeFunctionDerivsAtQuadPoints[jNode,:]
shapeFunctionQuadPointValuesJ = shapeFunctionAtQuadPoints[jNode,:]
quadraturePointValuesKinetic = (1.0/2.0)*shapeFunctionGradientQuadPointValuesI*shapeFunctionGradientQuadPointValuesJ
quadraturePointValuesPotX = (1/(normy*normz))*vx[start:end]*shapeFunctionQuadPointValuesI*shapeFunctionQuadPointValuesJ
quadraturePointValuesPotY = (1/(normx*normz))*vy[start:end]*shapeFunctionQuadPointValuesI*shapeFunctionQuadPointValuesJ
quadraturePointValuesPotZ = (1/(normx*normy))*vz[start:end]*shapeFunctionQuadPointValuesI*shapeFunctionQuadPointValuesJ
quadraturePointValuesMass = shapeFunctionQuadPointValuesI*shapeFunctionQuadPointValuesJ;
elementalInvJacAtQuadPoints = invJacobianQuadPointValues[start:end]
elementalJacAtQuadPoints = jacobianQuadPointValues[start:end]
shapeFunctionGradientIntegralIJ = fem.integrate(quadraturePointValuesKinetic,
elementalInvJacAtQuadPoints);
elementStiffnessMatX[iNode,jNode] = shapeFunctionGradientIntegralIJ +fem.integrate(quadraturePointValuesPotX,elementalJacAtQuadPoints)
elementStiffnessMatY[iNode,jNode] = shapeFunctionGradientIntegralIJ \
+fem.integrate(quadraturePointValuesPotY,elementalJacAtQuadPoints)
elementStiffnessMatZ[iNode,jNode] = shapeFunctionGradientIntegralIJ \
+fem.integrate(quadraturePointValuesPotZ,elementalJacAtQuadPoints)
elementMassMat[iNode,jNode] = fem.integrate(quadraturePointValuesMass,elementalJacAtQuadPoints);
#
# insert into global stiffness matrix
#
for iNode in range(numberNodesPerElement):
iGlobal = elementConnectivity[iNode,e]
for jNode in range(numberNodesPerElement):
jGlobal = elementConnectivity[jNode,e];
Hx[iGlobal,jGlobal] += elementStiffnessMatX[iNode,jNode]
Hy[iGlobal,jGlobal] += elementStiffnessMatY[iNode,jNode]
Hz[iGlobal,jGlobal] += elementStiffnessMatZ[iNode,jNode]
M[iGlobal,jGlobal] += elementMassMat[iNode,jNode];
Hx = self.condenseMatrix(Hx)
Hy = self.condenseMatrix(Hy)
Hz = self.condenseMatrix(Hz)
M = self.condenseMatrix(M)
return (Hx,Hy,Hz,M)