adaptation.py 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322
  1. import numpy as np
  2. from numpy import concatenate as cat
  3. import matplotlib.pyplot as plt
  4. import warnings
  5. from .preprocess import shape, discretization
  6. plt.rc('text', usetex=True)
  7. plt.rc('font', family='serif')
  8. # supress the deprecation warning
  9. warnings.filterwarnings("ignore", ".*GUI is implemented.*")
  10. class hdpg1d(object):
  11. """
  12. 1D HDG solver
  13. """
  14. def __init__(self, coeff):
  15. self.numEle = coeff.numEle
  16. self.numBasisFuncs = coeff.pOrder + 1
  17. self.tau_pos = coeff.tauPlus
  18. self.tau_neg = coeff.tauMinus
  19. self.c = coeff.convection
  20. self.kappa = coeff.diffusion
  21. self.coeff = coeff
  22. self.mesh = np.linspace(0, 1, self.numEle + 1)
  23. self.u = []
  24. self.estErrorList = [[], []]
  25. self.trueErrorList = [[], []]
  26. def bc(self, case, t=None):
  27. # boundary condition
  28. if case == 0:
  29. # advection-diffusion
  30. bc = [0, 0]
  31. if case == 1:
  32. # simple convection
  33. # bc = np.sin(2*np.pi*t)
  34. # adjoint boundary
  35. bc = [0, 1]
  36. return bc
  37. def forcing(self, x):
  38. # f = np.cos(2*np.pi*x)
  39. # f = 4*pi**2*sin(2*pi*x)
  40. f = 1
  41. return f
  42. def plotU(self, counter):
  43. """Plot solution u with smooth higher oredr quadrature"""
  44. uSmooth = np.array([])
  45. uNode = np.zeros(self.numEle + 1)
  46. xSmooth = np.array([])
  47. u = self.u[int(len(self.u) / 2):len(self.u)]
  48. # quadrature rule
  49. gorder = 10 * self.numBasisFuncs
  50. xi, wi = np.polynomial.legendre.leggauss(gorder)
  51. shp, shpx = shape(xi, self.numBasisFuncs)
  52. for j in range(1, self.numEle + 1):
  53. xSmooth = np.hstack((xSmooth, (self.mesh[(j - 1)] + self.mesh[j]) / 2 + (
  54. self.mesh[j] - self.mesh[j - 1]) / 2 * xi))
  55. uSmooth = np.hstack(
  56. (uSmooth, shp.T.dot(u[(j - 1) * self.numBasisFuncs:j * self.numBasisFuncs])))
  57. uNode[j - 1] = u[(j - 1) * self.numBasisFuncs]
  58. uNode[-1] = u[-1]
  59. plt.figure(1)
  60. plt.plot(xSmooth, uSmooth, '-', color='C3')
  61. plt.plot(self.mesh, uNode, 'C3.')
  62. plt.xlabel('$x$', fontsize=17)
  63. plt.ylabel('$u$', fontsize=17)
  64. # plt.axis([-0.05, 1.05, 0, 1.3])
  65. plt.grid()
  66. plt.pause(5e-1)
  67. plt.clf()
  68. def meshAdapt(self, index):
  69. """Given the index list, adapt the mesh"""
  70. inValue = np.zeros(len(index))
  71. for i in np.arange(len(index)):
  72. inValue[i] = (self.mesh[index[i]] +
  73. self.mesh[index[i] - 1]) / 2
  74. self.mesh = np.sort(np.insert(self.mesh, 0, inValue))
  75. def solveLocal(self):
  76. """Solve the primal problem"""
  77. if 'matLocal' in locals():
  78. # if matLocal exists,
  79. # only change the mesh instead of initializing again
  80. matLocal.mesh = self.mesh
  81. else:
  82. matLocal = discretization(self.coeff, self.mesh)
  83. F, L, R = matLocal.lhsGen()
  84. A, B, D = matLocal.eleGen()
  85. C, E, G, H = matLocal.interfaceGen()
  86. # solve
  87. K = -cat((C.T, G), axis=1)\
  88. .dot(np.linalg.inv(np.bmat([[A, -B], [B.T, D]]))
  89. .dot(cat((C, E)))) + H
  90. F_hat = np.array([L]).T - cat((C.T, G), axis=1)\
  91. .dot(np.linalg.inv(np.bmat([[A, -B], [B.T, D]])))\
  92. .dot(np.array([cat((R, F))]).T)
  93. uFace = np.linalg.solve(K, F_hat)
  94. u = np.linalg.inv(np.bmat([[A, -B], [B.T, D]]))\
  95. .dot(np.array([np.concatenate((R, F))]).T -
  96. cat((C, E)).dot(uFace))
  97. # self.u = u.A1
  98. return u.A1, uFace.A1
  99. def solveAdjoint(self):
  100. """Solve the adjoint problem"""
  101. # solve in the enriched space
  102. self.coeff.pOrder += 1
  103. if 'matAdjoint' in locals():
  104. matAdjoint.mesh = self.mesh
  105. else:
  106. matAdjoint = discretization(self.coeff, self.mesh)
  107. self.coeff.pOrder = self.coeff.pOrder - 1
  108. F, L, R = matAdjoint.lhsGen()
  109. A, B, D = matAdjoint.eleGen()
  110. C, E, G, H = matAdjoint.interfaceGen()
  111. # add adjoint LHS conditions
  112. F = np.zeros(len(F))
  113. R[-1] = -self.bc(1)[1]
  114. # assemble global matrix LHS
  115. LHS = np.bmat([[A, -B, C],
  116. [B.T, D, E],
  117. [C.T, G, H]])
  118. # solve
  119. U = np.linalg.solve(LHS.T, np.concatenate(
  120. (R, F, L)))
  121. return U[0:2 * len(C)], U[len(C):len(U)]
  122. def residual(self, U, hat_U, z, hat_z):
  123. if 'matResidual' in locals():
  124. # if matLocal exists,
  125. # only change the mesh instead of initializing again
  126. matResidual.mesh = self.mesh
  127. else:
  128. matResidual = discretization(self.coeff, self.mesh, 1)
  129. F, L, R = matResidual.lhsGen()
  130. # A, B, D = matResidual.eleGen()
  131. # C, E, G, H = matResidual.interfaceGen()
  132. numEle = self.numEle
  133. p = self.numBasisFuncs + 1
  134. p_l = p - 1
  135. # order of gauss quadrature
  136. gorder = 2 * p
  137. # shape function and gauss quadrature
  138. xi, wi = np.polynomial.legendre.leggauss(gorder)
  139. shp, shpx = shape(xi, p)
  140. shp_l, shpx_l = shape(xi, p_l)
  141. # ---------------------------------------------------------------------
  142. # advection constant
  143. con = self.c
  144. # diffusion constant
  145. kappa = self.kappa
  146. z_q, z_u, z_hat = np.zeros(p * numEle), \
  147. np.zeros(p *
  148. numEle), np.zeros(numEle - 1)
  149. q, u, lamba = np.zeros(p_l * numEle), \
  150. np.zeros(p_l * numEle), np.zeros(numEle - 1)
  151. for i in np.arange(p * numEle):
  152. z_q[i] = z[i]
  153. z_u[i] = z[i + p * numEle]
  154. for i in np.arange(p_l * numEle):
  155. q[i] = U[i]
  156. u[i] = U[i + p_l * numEle]
  157. for i in np.arange(numEle - 1):
  158. z_hat[i] = hat_z[i]
  159. # add boundary condtions to U_hat
  160. U_hat = np.zeros(numEle + 1)
  161. for i, x in enumerate(hat_U):
  162. U_hat[i + 1] = x
  163. U_hat[0] = self.bc(0)[0]
  164. U_hat[-1] = self.bc(0)[1]
  165. # L, easy in 1d
  166. # L = np.zeros(numEle + 1)
  167. # R, easy in 1d
  168. RR = np.zeros(p * numEle)
  169. # elemental forcing vector
  170. dist = np.zeros(numEle)
  171. for i in range(1, numEle + 1):
  172. dist[i - 1] = self.mesh[i] - self.mesh[i - 1]
  173. # elemental h
  174. h = np.zeros((2, 2))
  175. h[0, 0], h[-1, -1] = -con - self.tau_pos, con - self.tau_neg
  176. # mappinng matrix
  177. map_h = np.zeros((2, numEle), dtype=int)
  178. map_h[:, 0] = np.arange(2)
  179. for i in np.arange(1, numEle):
  180. map_h[:, i] = np.arange(
  181. map_h[2 - 1, i - 1], map_h[2 - 1, i - 1] + 2)
  182. # assemble H and eliminate boundaries
  183. H = np.zeros((numEle + 1, numEle + 1))
  184. for i in range(numEle):
  185. for j in range(2):
  186. m = map_h[j, i]
  187. for k in range(2):
  188. n = map_h[k, i]
  189. H[m, n] += h[j, k]
  190. H = H[1:numEle][:, 1:numEle]
  191. # elemental g
  192. g = np.zeros((2, p_l))
  193. g[0, 0], g[-1, -1] = self.tau_pos, self.tau_neg
  194. # mapping matrix
  195. map_g_x = map_h
  196. map_g_y = np.arange(p_l * numEle, dtype=int).reshape(numEle, p_l).T
  197. # assemble global G
  198. G = np.zeros((numEle + 1, p_l * numEle))
  199. for i in range(numEle):
  200. for j in range(2):
  201. m = map_g_x[j, i]
  202. for k in range(p_l):
  203. n = map_g_y[k, i]
  204. G[m, n] += g[j, k]
  205. G = G[1:numEle, :]
  206. # elemental c
  207. c = np.zeros((p_l, 2))
  208. c[0, 0], c[-1, -1] = -1, 1
  209. # mapping matrix
  210. map_e_x = np.arange(p_l * numEle, dtype=int).reshape(numEle, p_l).T
  211. map_e_y = map_h
  212. # assemble global C
  213. C = np.zeros((p_l * numEle, numEle + 1))
  214. for i in range(numEle):
  215. for j in range(p_l):
  216. m = map_e_x[j, i]
  217. for k in range(2):
  218. n = map_e_y[k, i]
  219. C[m, n] += c[j, k]
  220. C = C[:, 1:numEle]
  221. # L, easy in 1d
  222. # L = np.zeros(numEle - 1)
  223. # residual vector
  224. Residual = np.zeros(self.numEle)
  225. for i in np.arange(self.numEle):
  226. a = dist[i] / 2 * 1 / kappa * \
  227. ((shp.T).T).dot(np.diag(wi).dot(shp_l.T))
  228. b = ((shpx.T) * np.ones((gorder, p))
  229. ).T.dot(np.diag(wi).dot(shp_l.T))
  230. b_t = ((shpx_l.T) * np.ones((gorder, p_l))
  231. ).T.dot(np.diag(wi).dot(shp.T))
  232. d = dist[i] / 2 * shp.dot(np.diag(wi).dot(shp_l.T))
  233. d[0, 0] += self.tau_pos
  234. d[-1, -1] += self.tau_neg
  235. h = np.zeros((2, 2))
  236. h[0, 0], h[-1, -1] = -con - self.tau_pos, con - self.tau_neg
  237. g = np.zeros((2, p_l))
  238. g[0, 0], g[-1, -1] = self.tau_pos, self.tau_neg
  239. e = np.zeros((p, 2))
  240. e[0, 0], e[-1, -1] = -con - self.tau_pos, con - self.tau_neg
  241. c = np.zeros((p, 2))
  242. c[0, 0], c[-1, -1] = -1, 1
  243. m = np.zeros((2, p_l))
  244. m[0, 0], m[-1, -1] = -1, 1
  245. # local error
  246. Residual[i] = (np.concatenate((a.dot(q[p_l * i:p_l * i + p_l]) + -b.dot(u[p_l * i:p_l * i + p_l]) + c.dot(U_hat[i:i + 2]),
  247. b_t.T.dot(q[p_l * i:p_l * i + p_l]) + d.dot(u[p_l * i:p_l * i + p_l]) + e.dot(U_hat[i:i + 2]))) - np.concatenate((R[p * i:p * i + p], F[p * i:p * i + p]))).dot(1 - np.concatenate((z_q[p * i:p * i + p], z_u[p * i:p * i + p])))
  248. com_index = np.argsort(np.abs(Residual))
  249. # select \theta% elements with the large error
  250. theta = 0.15
  251. refine_index = com_index[int(self.numEle * (1 - theta)):len(Residual)]
  252. # global error
  253. R_g = (C.T.dot(q) + G.dot(u) + H.dot(U_hat[1:-1])).dot(1 - z_hat)
  254. return np.abs(np.sum(Residual) + np.sum(R_g)), refine_index + 1
  255. def adaptive(self):
  256. tol = 1e-10
  257. estError = 10
  258. counter = 0
  259. ceilCounter = 30
  260. while estError > tol and counter < ceilCounter:
  261. print("Iteration {}. Target function error {:.3e}.".format(
  262. counter, estError))
  263. # solve
  264. u, uFace = self.solveLocal()
  265. adjoint, adjointFace = self.solveAdjoint()
  266. self.u = u
  267. # plot the solution at certain counter
  268. if counter in [0, 4, 9, 19]:
  269. self.plotU(counter)
  270. # record error
  271. self.trueErrorList[0].append(self.numEle)
  272. self.trueErrorList[1].append(
  273. u[self.numEle * self.numBasisFuncs - 1])
  274. estError, index = self.residual(u, uFace, adjoint, adjointFace)
  275. self.estErrorList[0].append(self.numEle)
  276. self.estErrorList[1].append(estError)
  277. # adapt
  278. index = index.tolist()
  279. self.meshAdapt(index)
  280. self.numEle = self.numEle + len(index)
  281. counter += 1