preprocess.py 8.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227
  1. import json
  2. import numpy as np
  3. import os
  4. from collections import namedtuple
  5. from scipy.linalg import block_diag
  6. for loc in os.curdir, os.path.expanduser("~"), "/etc/hdpg1d":
  7. try:
  8. with open(os.path.join(loc, "config.json")) as source:
  9. data = json.load(source)
  10. except IOError:
  11. pass
  12. def shape(x, p):
  13. """generate p shape functions and its first order derivative
  14. at the given location x. x can be an array"""
  15. A = np.array([np.linspace(-1, 1, p)]).T**np.arange(p)
  16. C = np.linalg.inv(A).T
  17. x = np.array([x]).T
  18. shp = C.dot((x**np.arange(p)).T)
  19. shpx = C[:, 1::1].dot((x**np.arange(p - 1) * np.arange(1, p)).T)
  20. return shp, shpx
  21. def forcing(x):
  22. f = 1
  23. return f
  24. def boundaryCondition(case):
  25. if case == 0:
  26. # primal problem
  27. bc = [0, 0]
  28. if case == 1:
  29. # adjoint problem
  30. bc = [0, 1]
  31. return bc
  32. class discretization(object):
  33. """Given the problem statement and current mesh,
  34. construct the discretization matrices"""
  35. def __init__(self, coeff, mesh, enrich=None):
  36. self.mesh = mesh
  37. self.coeff = coeff
  38. self.enrich = enrich
  39. # the following init are for the sake of simplicity
  40. self.numBasisFuncs = coeff.pOrder + 1
  41. self.tau_pos = coeff.TAUPLUS
  42. self.tau_neg = coeff.TAUMINUS
  43. self.conv = coeff.CONVECTION
  44. self.kappa = coeff.DIFFUSION
  45. self.numEle = len(mesh) - 1
  46. self.dist = self.distGen()
  47. # shape function and gauss quadrature
  48. self.gqOrder = 5 * self.numBasisFuncs
  49. self.xi, self.wi = np.polynomial.legendre.leggauss(self.gqOrder)
  50. self.shp, self.shpx = shape(self.xi, self.numBasisFuncs)
  51. # enrich the space if the enrich argument is given
  52. if enrich is not None:
  53. self.shpEnrich, self.shpxEnrich = shape(
  54. self.xi, self.numBasisFuncs + enrich)
  55. def distGen(self):
  56. dist = np.zeros(self.numEle)
  57. for i in range(1, self.numEle + 1):
  58. dist[i - 1] = self.mesh[i] - self.mesh[i - 1]
  59. return dist
  60. def matGroup(self):
  61. lhsMat = self.lhsGen()
  62. F, L, R = lhsMat.F, lhsMat.L, lhsMat.R
  63. eleMat = self.eleGen()
  64. A, B, BonU, D = eleMat.A, eleMat.B, eleMat.BonU, eleMat.D
  65. faceMat = self.interfaceGen()
  66. C, E, G, H = faceMat.C, faceMat.E, faceMat.G, faceMat.H
  67. matGroup = namedtuple(
  68. 'matGroup', ['A', 'B', 'BonU', 'C', 'D', 'E', 'F', 'G', 'H', 'L', 'R'])
  69. return matGroup(A, B, BonU, C, D, E, F, G, H, L, R)
  70. def lhsGen(self):
  71. """Generate matrices associated with left hand side"""
  72. if self.enrich is not None:
  73. # when enrich is given
  74. # provide matrices used in the DWR residual calculation
  75. numBasisFuncs = self.numBasisFuncs + self.enrich
  76. shp = self.shpEnrich
  77. else:
  78. numBasisFuncs = self.numBasisFuncs
  79. shp = self.shp
  80. # forcing vector F
  81. F = np.zeros(numBasisFuncs * self.numEle)
  82. for i in range(1, self.numEle + 1):
  83. f = self.dist[i - 1] / 2 \
  84. * shp.dot(self.wi * forcing(self.mesh[i - 1] +
  85. 1 / 2 * (1 + self.xi) *
  86. self.dist[i - 1]))
  87. F[(i - 1) * numBasisFuncs:i * numBasisFuncs] = f
  88. F[0] += (self.conv + self.tau_pos) * boundaryCondition(0)[0]
  89. F[-1] += (-self.conv + self.tau_neg) * boundaryCondition(0)[1]
  90. # L, easy in 1d
  91. L = np.zeros(self.numEle - 1)
  92. # R, easy in 1d
  93. R = np.zeros(numBasisFuncs * self.numEle)
  94. R[0] = boundaryCondition(0)[0]
  95. R[-1] = -boundaryCondition(0)[1]
  96. lhsMat = namedtuple('lhsMat', ['F', 'L', 'R'])
  97. return lhsMat(F, L, R)
  98. def eleGen(self):
  99. """Generate matrices associated with elements"""
  100. if self.enrich is not None:
  101. numBasisFuncsEnrich = self.numBasisFuncs + self.enrich
  102. shpEnrich = self.shpEnrich
  103. shpxEnrich = self.shpxEnrich
  104. b = (self.shpx.T * np.ones((self.gqOrder, self.numBasisFuncs))
  105. ).T.dot(np.diag(self.wi).dot(shpEnrich.T))
  106. # BonQ is only used in calculating DWR residual
  107. BonQ = block_diag(*[b] * (self.numEle)).T
  108. else:
  109. numBasisFuncsEnrich = self.numBasisFuncs
  110. shpEnrich = self.shp
  111. shpxEnrich = self.shpx
  112. BonQ = None
  113. a = 1 / self.coeff.DIFFUSION * \
  114. shpEnrich.dot(np.diag(self.wi).dot(self.shp.T))
  115. A = np.repeat(self.dist, self.numBasisFuncs) / \
  116. 2 * block_diag(*[a] * (self.numEle))
  117. b = (shpxEnrich.T * np.ones((self.gqOrder, numBasisFuncsEnrich))
  118. ).T.dot(np.diag(self.wi).dot(self.shp.T))
  119. B = block_diag(*[b] * (self.numEle))
  120. d = self.coeff.REACTION * \
  121. shpEnrich.dot(np.diag(self.wi).dot(self.shp.T))
  122. # assemble D
  123. dFace = np.zeros((numBasisFuncsEnrich, self.numBasisFuncs))
  124. dFace[0, 0] = self.tau_pos
  125. dFace[-1, -1] = self.tau_neg
  126. dConv = -self.conv * (shpxEnrich.T * np.ones((self.gqOrder,
  127. numBasisFuncsEnrich)))\
  128. .T.dot(np.diag(self.wi).dot(self.shp.T))
  129. D = np.repeat(self.dist, self.numBasisFuncs) / 2 * \
  130. block_diag(*[d] * (self.numEle)) + \
  131. block_diag(*[dFace] * (self.numEle)) +\
  132. block_diag(*[dConv] * (self.numEle))
  133. eleMat = namedtuple('eleMat', ['A', 'B', 'BonU', 'D'])
  134. return eleMat(A, B, BonQ, D)
  135. def interfaceGen(self):
  136. """Generate matrices associated with interfaces"""
  137. if self.enrich is not None:
  138. # when enrich is given
  139. # provide matrices used in the DWR residual calculation
  140. numBasisFuncs = self.numBasisFuncs + self.enrich
  141. else:
  142. numBasisFuncs = self.numBasisFuncs
  143. tau_pos, tau_neg, conv = self.tau_pos, self.tau_neg, self.conv
  144. # elemental h
  145. h = np.zeros((2, 2))
  146. h[0, 0], h[-1, -1] = -conv - tau_pos, conv - tau_neg
  147. # mappinng matrix
  148. map_h = np.zeros((2, self.numEle), dtype=int)
  149. map_h[:, 0] = np.arange(2)
  150. for i in np.arange(1, self.numEle):
  151. map_h[:, i] = np.arange(
  152. map_h[2 - 1, i - 1], map_h[2 - 1, i - 1] + 2)
  153. # assemble H and eliminate boundaries
  154. H = np.zeros((self.numEle + 1, self.numEle + 1))
  155. for i in range(self.numEle):
  156. for j in range(2):
  157. m = map_h[j, i]
  158. for k in range(2):
  159. n = map_h[k, i]
  160. H[m, n] += h[j, k]
  161. H = H[1:self.numEle][:, 1:self.numEle]
  162. # elemental e
  163. e = np.zeros((numBasisFuncs, 2))
  164. e[0, 0], e[-1, -1] = -conv - tau_pos, conv - tau_neg
  165. # mapping matrix
  166. map_e_x = np.arange(numBasisFuncs * self.numEle,
  167. dtype=int).reshape(self.numEle,
  168. numBasisFuncs).T
  169. map_e_y = map_h
  170. # assemble global E
  171. E = np.zeros((numBasisFuncs * self.numEle, self.numEle + 1))
  172. for i in range(self.numEle):
  173. for j in range(numBasisFuncs):
  174. m = map_e_x[j, i]
  175. for k in range(2):
  176. n = map_e_y[k, i]
  177. E[m, n] += e[j, k]
  178. E = E[:, 1:self.numEle]
  179. # elemental c
  180. c = np.zeros((numBasisFuncs, 2))
  181. c[0, 0], c[-1, -1] = -1, 1
  182. # assemble global C
  183. C = np.zeros((numBasisFuncs * self.numEle, self.numEle + 1))
  184. for i in range(self.numEle):
  185. for j in range(numBasisFuncs):
  186. m = map_e_x[j, i]
  187. for k in range(2):
  188. n = map_e_y[k, i]
  189. C[m, n] += c[j, k]
  190. C = C[:, 1:self.numEle]
  191. # elemental g
  192. g = np.zeros((2, numBasisFuncs))
  193. g[0, 0], g[-1, -1] = tau_pos, tau_neg
  194. # mapping matrix
  195. map_g_x = map_h
  196. map_g_y = np.arange(numBasisFuncs * self.numEle,
  197. dtype=int).reshape(self.numEle,
  198. numBasisFuncs).T
  199. # assemble global G
  200. G = np.zeros((self.numEle + 1, numBasisFuncs * self.numEle))
  201. for i in range(self.numEle):
  202. for j in range(2):
  203. m = map_g_x[j, i]
  204. for k in range(numBasisFuncs):
  205. n = map_g_y[k, i]
  206. G[m, n] += g[j, k]
  207. G = G[1:self.numEle, :]
  208. faceMat = namedtuple('faceMat', ['C', 'E', 'G', 'H'])
  209. return faceMat(C, E, G, H)