discretization.py 19 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554
  1. import numpy as np
  2. from numpy import sin, cos, pi
  3. from scipy.linalg import block_diag
  4. import matplotlib.pyplot as plt
  5. class HDPG1d(object):
  6. """
  7. 1d advection hdg solver outlined in 'an implicit HHDG method for
  8. confusion'. Test case: /tau = 1, convection only, linear and higher order.
  9. Please enter number of elements and polynomial order, i.e., HDG1d(10,2)
  10. """
  11. # stablization parameter
  12. tau_pos = 1e-6
  13. tau_neg = 1e-6
  14. # advection constant
  15. c = 0
  16. # diffusion constant
  17. kappa = 1e-6
  18. def __init__(self, n_ele, p):
  19. self.n_ele = n_ele
  20. self.p = p + 1
  21. self.estError = []
  22. self.trueError = []
  23. def bc(self, case, t=None):
  24. # boundary condition
  25. if case == 0:
  26. # advection-diffusion
  27. bc = [0, 0]
  28. if case == 1:
  29. # simple convection
  30. # bc = np.sin(2*np.pi*t)
  31. # adjoint boundary
  32. bc = [0, 1]
  33. return bc
  34. def shape(self, x, p):
  35. """ evaluate shape functions at give locations"""
  36. # coeffient matrix
  37. A = np.array([np.linspace(-1, 1, p)]).T**np.arange(p)
  38. C = np.linalg.inv(A).T
  39. x = np.array([x]).T
  40. shp = C.dot((x**np.arange(p)).T)
  41. shpx = C[:, 1::1].dot((x**np.arange(p - 1) * np.arange(1, p)).T)
  42. return shp, shpx
  43. def forcing(self, x):
  44. # f = np.cos(2*np.pi*x)
  45. # f = 4*pi**2*sin(2*pi*x)
  46. f = 1
  47. return f
  48. def mesh(self, n_ele, index, x):
  49. """generate mesh"""
  50. # if n_ele < 1 or n_ele > self.n_ele:
  51. # raise RuntimeError('Bad Element number')
  52. in_value = np.zeros(len(index))
  53. for i in np.arange(len(index)):
  54. in_value[i] = (x[index[i]] + x[index[i] - 1]) / 2
  55. x_c = np.sort(np.insert(x, 0, in_value))
  56. x_i = np.linspace(x_c[n_ele - 1], x_c[n_ele], num=self.p)
  57. dx = x_c[n_ele] - x_c[n_ele - 1]
  58. return x_i, dx, x_c
  59. def uexact(self, x):
  60. # Ue = self.bc(0) + np.sin(2*pi*x)
  61. Ue = (np.exp(1 / self.kappa * (x - 1)) - 1) / \
  62. (np.exp(-1 / self.kappa) - 1)
  63. return Ue
  64. def matrix_gen(self, index, x):
  65. n_ele = self.n_ele
  66. # # space size
  67. # dx = 1/n_ele
  68. # order of polynomial shape functions
  69. p = self.p
  70. # order of gauss quadrature
  71. gorder = 2 * p
  72. # shape function and gauss quadrature
  73. xi, wi = np.polynomial.legendre.leggauss(gorder)
  74. shp, shpx = self.shape(xi, p)
  75. # ---------------------------------------------------------------------
  76. # advection constant
  77. con = self.c
  78. # diffusion constant
  79. kappa = self.kappa
  80. # number of nodes (solution U)
  81. n_ele = self.n_ele + len(index)
  82. # elemental forcing vector
  83. F = np.zeros(p * n_ele)
  84. for i in range(1, n_ele + 1):
  85. x_i, dx_i, _ = self.mesh(i, index, x)
  86. f = dx_i / 2 * \
  87. shp.dot(wi * self.forcing(x[0] + 1 / 2 * (1 + xi) * dx_i))
  88. F[(i - 1) * p:(i - 1) * p + p] = f
  89. # elemental m (for convection only)
  90. # m = dx/2*shp.dot(np.diag(wi).dot(shp.T))
  91. # add boundary
  92. F[0] += (con + self.tau_pos) * self.bc(0)[0]
  93. F[-1] += (-con + self.tau_neg) * self.bc(0)[1]
  94. # elemental d
  95. # d = con*(-shpx.T*np.ones((gorder,p))).T.dot(np.diag(wi).dot(shp.T))
  96. d = shp.dot(np.diag(wi).dot(shp.T))
  97. # d[0, 0] += self.tau_pos
  98. # d[-1, -1] += self.tau_neg
  99. # elemental a
  100. a = 1 / kappa * shp.dot(np.diag(wi).dot(shp.T))
  101. # elemental b
  102. b = (shpx.T * np.ones((gorder, p))).T.dot(np.diag(wi).dot(shp.T))
  103. # elemental h
  104. h = np.zeros((2, 2))
  105. h[0, 0], h[-1, -1] = -con - self.tau_pos, con - self.tau_neg
  106. # mappinng matrix
  107. map_h = np.zeros((2, n_ele), dtype=int)
  108. map_h[:, 0] = np.arange(2)
  109. for i in np.arange(1, n_ele):
  110. map_h[:, i] = np.arange(
  111. map_h[2 - 1, i - 1], map_h[2 - 1, i - 1] + 2)
  112. # assemble H and eliminate boundaries
  113. H = np.zeros((n_ele + 1, n_ele + 1))
  114. for i in range(n_ele):
  115. for j in range(2):
  116. m = map_h[j, i]
  117. for k in range(2):
  118. n = map_h[k, i]
  119. H[m, n] += h[j, k]
  120. H = H[1:n_ele][:, 1:n_ele]
  121. # elemental g
  122. g = np.zeros((2, p))
  123. g[0, 0], g[-1, -1] = self.tau_pos, self.tau_neg
  124. # mapping matrix
  125. map_g_x = map_h
  126. map_g_y = np.arange(p * n_ele, dtype=int).reshape(n_ele, p).T
  127. # assemble global G
  128. G = np.zeros((n_ele + 1, p * n_ele))
  129. for i in range(n_ele):
  130. for j in range(2):
  131. m = map_g_x[j, i]
  132. for k in range(p):
  133. n = map_g_y[k, i]
  134. G[m, n] += g[j, k]
  135. G = G[1:n_ele, :]
  136. # elemental e
  137. e = np.zeros((p, 2))
  138. e[0, 0], e[-1, -1] = -con - self.tau_pos, con - self.tau_neg
  139. # mapping matrix
  140. map_e_x = np.arange(p * n_ele, dtype=int).reshape(n_ele, p).T
  141. map_e_y = map_h
  142. # assemble global E
  143. E = np.zeros((p * n_ele, n_ele + 1))
  144. for i in range(n_ele):
  145. for j in range(p):
  146. m = map_e_x[j, i]
  147. for k in range(2):
  148. n = map_e_y[k, i]
  149. E[m, n] += e[j, k]
  150. E = E[:, 1:n_ele]
  151. # L, easy in 1d
  152. L = np.zeros(n_ele - 1)
  153. # elemental c
  154. c = np.zeros((p, 2))
  155. c[0, 0], c[-1, -1] = -1, 1
  156. # assemble global C
  157. C = np.zeros((p * n_ele, n_ele + 1))
  158. for i in range(n_ele):
  159. for j in range(p):
  160. m = map_e_x[j, i]
  161. for k in range(2):
  162. n = map_e_y[k, i]
  163. C[m, n] += c[j, k]
  164. C = C[:, 1:n_ele]
  165. # L, easy in 1d
  166. L = np.zeros(n_ele - 1)
  167. # R, easy in 1d
  168. R = np.zeros(p * n_ele)
  169. R[0] = self.bc(0)[0]
  170. R[-1] = -self.bc(0)[1]
  171. return d, m, E, G, H, F, L, a, b, C, R
  172. def solve_local(self, index, x):
  173. """ solve the 1d advection equation wit local HDG"""
  174. d, _, E, G, H, F, L, a, b, C, R = self.matrix_gen(index, x)
  175. # find dx
  176. dx = np.zeros(self.n_ele + len(index))
  177. for i in range(1, self.n_ele + len(index) + 1):
  178. x_i, dx_i, x_n = self.mesh(i, index, x)
  179. dx[i - 1] = dx_i
  180. # assemble global D
  181. bb = np.zeros((self.p, self.p))
  182. bb[0, 0] = self.tau_pos
  183. bb[-1, -1] = self.tau_neg
  184. D = np.repeat(dx, self.p) / 2 * block_diag(*[d] * (
  185. self.n_ele + len(index))) + block_diag(*[bb] * (self.n_ele + len(index)))
  186. # assemble global A
  187. A = np.repeat(dx, self.p) / 2 * block_diag(*
  188. [a] * (self.n_ele + len(index)))
  189. # assemble global B
  190. B = block_diag(*[b] * (self.n_ele + len(index)))
  191. # solve U and \lambda
  192. K = -np.concatenate((C.T, G), axis=1).dot(np.linalg.inv(
  193. np.bmat([[A, -B], [B.T, D]])).dot(np.concatenate((C, E)))) + H
  194. F_hat = np.array([L]).T - np.concatenate((C.T, G), axis=1).dot(np.linalg.inv(
  195. np.bmat([[A, -B], [B.T, D]]))).dot(np.array([np.concatenate((R, F))]).T)
  196. lamba = np.linalg.solve(K, F_hat)
  197. U = np.linalg.inv(np.bmat([[A, -B], [B.T, D]])).dot(
  198. np.array([np.concatenate((R, F))]).T - np.concatenate((C, E)).dot(lamba))
  199. return U, lamba
  200. def solve_adjoint(self, index, x, u, u_hat):
  201. self.p = self.p + 1
  202. d, _, E, G, H, F, L, a, b, C, R = self.matrix_gen(index, x)
  203. # add boundary
  204. F = np.zeros(len(F))
  205. R[-1] = -self.bc(1)[1]
  206. # find dx
  207. dx = np.zeros(self.n_ele + len(index))
  208. for i in range(1, self.n_ele + len(index) + 1):
  209. x_i, dx_i, x_n = self.mesh(i, index, x)
  210. dx[i - 1] = dx_i
  211. # assemble global D
  212. bb = np.zeros((self.p, self.p))
  213. bb[0, 0] = self.tau_pos
  214. bb[-1, -1] = self.tau_neg
  215. D = np.repeat(dx, self.p) / 2 * block_diag(*[d] * (
  216. self.n_ele + len(index))) + block_diag(*[bb] * (self.n_ele + len(index)))
  217. # assemble global A
  218. A = np.repeat(dx, self.p) / 2 * block_diag(*
  219. [a] * (self.n_ele + len(index)))
  220. # assemble global B
  221. B = block_diag(*[b] * (self.n_ele + len(index)))
  222. # # assemble global matrix LHS
  223. LHS = np.bmat([[A, -B, C], [B.T, D, E], [C.T, G, H]])
  224. # solve U and \lambda
  225. U = np.linalg.solve(LHS.T, np.concatenate((R, F, L)))
  226. return U[0:2 * self.p * (self.n_ele + len(index))], U[2 * self.p * (self.n_ele + len(index)):len(U)]
  227. def diffusion(self):
  228. """solve 1d convection with local HDG"""
  229. # begin and end time
  230. t, T = 0, 1
  231. # time marching step for diffusion equation
  232. dt = 1e-3
  233. d, m, E, G, H, F, L, a, b, C, R = self.matrix_gen()
  234. # add time derivatives to the space derivatives (both are
  235. # elmental-wise)
  236. d = d + 1 / dt * m
  237. # assemble global D
  238. D = block_diag(*[d] * self.n_ele)
  239. # assemble global A
  240. A = block_diag(*[a] * self.n_ele)
  241. # assemble global B
  242. B = block_diag(*[b] * self.n_ele)
  243. # initial condition
  244. X = np.zeros(self.p * self.n_ele)
  245. for i in range(1, self.n_ele + 1):
  246. x = self.mesh(i)
  247. X[(i - 1) * self.p:(i - 1) * self.p + self.p] = x
  248. U = np.concatenate((pi * cos(pi * X), sin(pi * X)))
  249. # assemble M
  250. M = block_diag(*[1 / dt * m] * self.n_ele)
  251. # time marching
  252. while t < T:
  253. # add boundary conditions
  254. F_dynamic = F + \
  255. M.dot(U[self.n_ele * self.p:2 * self.n_ele * self.p])
  256. # assemble global matrix LHS
  257. LHS = np.bmat([[A, -B, C], [B.T, D, E], [C.T, G, H]])
  258. # solve U and \lambda
  259. U = np.linalg.solve(LHS, np.concatenate((R, F_dynamic, L)))
  260. # plot solutions
  261. plt.clf()
  262. plt.plot(X, U[self.n_ele * self.p:2 * self.n_ele * self.p], '-r.')
  263. plt.plot(X, sin(pi * X) * np.exp(-pi**2 * t))
  264. plt.ylim([0, 1])
  265. plt.grid()
  266. plt.pause(1e-3)
  267. plt.close()
  268. print("Diffusion equation du/dt - du^2/d^2x = 0 with u_exact ="
  269. ' 6sin(pi*x)*exp(-pi^2*t).')
  270. plt.figure(1)
  271. plt.plot(X, U[self.n_ele * self.p:2 * self.n_ele * self.p], '-r.')
  272. plt.plot(X, sin(pi * X) * np.exp(-pi**2 * T))
  273. plt.xlabel('x')
  274. plt.ylabel('y')
  275. plt.legend(('Numberical', 'Exact'), loc='upper left')
  276. plt.title('Simple Diffusion Equation Solution at t = {}'.format(T))
  277. plt.grid()
  278. plt.savefig('diffusion', bbox_inches='tight')
  279. plt.show(block=False)
  280. return U
  281. def residual(self, U, hat_U, z, hat_z, dx, index, x_c):
  282. n_ele = self.n_ele
  283. # order of polynomial shape functions
  284. p = self.p
  285. p_l = p - 1
  286. # order of gauss quadrature
  287. gorder = 2 * p
  288. # shape function and gauss quadrature
  289. xi, wi = np.polynomial.legendre.leggauss(gorder)
  290. shp, shpx = self.shape(xi, p)
  291. shp_l, shpx_l = self.shape(xi, p_l)
  292. # ---------------------------------------------------------------------
  293. # advection constant
  294. con = self.c
  295. # diffusion constant
  296. kappa = self.kappa
  297. z_q, z_u, z_hat = np.zeros(self.p * self.n_ele), \
  298. np.zeros(self.p * self.n_ele), np.zeros(self.n_ele - 1)
  299. q, u, lamba = np.zeros(p_l * self.n_ele), \
  300. np.zeros(p_l * self.n_ele), np.zeros(self.n_ele - 1)
  301. for i in np.arange(self.p * self.n_ele):
  302. z_q[i] = z[i]
  303. z_u[i] = z[i + self.p * self.n_ele]
  304. for i in np.arange(p_l * self.n_ele):
  305. q[i] = U[i]
  306. u[i] = U[i + p_l * self.n_ele]
  307. for i in np.arange(self.n_ele - 1):
  308. z_hat[i] = hat_z[i]
  309. # add boundary condtions to U_hat
  310. U_hat = np.zeros(self.n_ele + 1)
  311. for i, x in enumerate(hat_U):
  312. U_hat[i + 1] = x
  313. U_hat[0] = self.bc(0)[0]
  314. U_hat[-1] = self.bc(0)[1]
  315. # L, easy in 1d
  316. L = np.zeros(n_ele + 1)
  317. # R, easy in 1d
  318. RR = np.zeros(p * n_ele)
  319. # elemental forcing vector
  320. F = np.zeros(p * n_ele)
  321. for i in range(1, n_ele + 1):
  322. f = dx[i - 1] / 2 * \
  323. shp.dot(
  324. wi * self.forcing(x_c[0] + 1 / 2 * (1 + xi) * dx[i - 1]))
  325. F[(i - 1) * p:(i - 1) * p + p] = f
  326. # elemental h
  327. h = np.zeros((2, 2))
  328. h[0, 0], h[-1, -1] = -con - self.tau_pos, con - self.tau_neg
  329. # mappinng matrix
  330. map_h = np.zeros((2, n_ele), dtype=int)
  331. map_h[:, 0] = np.arange(2)
  332. for i in np.arange(1, n_ele):
  333. map_h[:, i] = np.arange(
  334. map_h[2 - 1, i - 1], map_h[2 - 1, i - 1] + 2)
  335. # assemble H and eliminate boundaries
  336. H = np.zeros((n_ele + 1, n_ele + 1))
  337. for i in range(n_ele):
  338. for j in range(2):
  339. m = map_h[j, i]
  340. for k in range(2):
  341. n = map_h[k, i]
  342. H[m, n] += h[j, k]
  343. H = H[1:n_ele][:, 1:n_ele]
  344. # elemental g
  345. g = np.zeros((2, p_l))
  346. g[0, 0], g[-1, -1] = self.tau_pos, self.tau_neg
  347. # mapping matrix
  348. map_g_x = map_h
  349. map_g_y = np.arange(p_l * n_ele, dtype=int).reshape(n_ele, p_l).T
  350. # assemble global G
  351. G = np.zeros((n_ele + 1, p_l * n_ele))
  352. for i in range(n_ele):
  353. for j in range(2):
  354. m = map_g_x[j, i]
  355. for k in range(p_l):
  356. n = map_g_y[k, i]
  357. G[m, n] += g[j, k]
  358. G = G[1:n_ele, :]
  359. # elemental c
  360. c = np.zeros((p_l, 2))
  361. c[0, 0], c[-1, -1] = -1, 1
  362. # mapping matrix
  363. map_e_x = np.arange(p_l * n_ele, dtype=int).reshape(n_ele, p_l).T
  364. map_e_y = map_h
  365. # assemble global C
  366. C = np.zeros((p_l * n_ele, n_ele + 1))
  367. for i in range(n_ele):
  368. for j in range(p_l):
  369. m = map_e_x[j, i]
  370. for k in range(2):
  371. n = map_e_y[k, i]
  372. C[m, n] += c[j, k]
  373. C = C[:, 1:n_ele]
  374. # L, easy in 1d
  375. L = np.zeros(n_ele - 1)
  376. # residual vector
  377. R = np.zeros(self.n_ele)
  378. for i in np.arange(self.n_ele):
  379. a = dx[i] / 2 * 1 / kappa * \
  380. ((shp.T).T).dot(np.diag(wi).dot(shp_l.T))
  381. b = ((shpx.T) * np.ones((gorder, p))
  382. ).T.dot(np.diag(wi).dot(shp_l.T))
  383. b_t = ((shpx_l.T) * np.ones((gorder, p_l))
  384. ).T.dot(np.diag(wi).dot(shp.T))
  385. d = dx[i] / 2 * shp.dot(np.diag(wi).dot(shp_l.T))
  386. d[0, 0] += self.tau_pos
  387. d[-1, -1] += self.tau_neg
  388. h = np.zeros((2, 2))
  389. h[0, 0], h[-1, -1] = -con - self.tau_pos, con - self.tau_neg
  390. g = np.zeros((2, p_l))
  391. g[0, 0], g[-1, -1] = self.tau_pos, self.tau_neg
  392. e = np.zeros((p, 2))
  393. e[0, 0], e[-1, -1] = -con - self.tau_pos, con - self.tau_neg
  394. c = np.zeros((p, 2))
  395. c[0, 0], c[-1, -1] = -1, 1
  396. m = np.zeros((2, p_l))
  397. m[0, 0], m[-1, -1] = -1, 1
  398. # local error
  399. R[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]),
  400. 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((RR[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])))
  401. com_index = np.argsort(np.abs(R))
  402. # select \theta% elements with the large error
  403. theta = 0.15
  404. refine_index = com_index[int(self.n_ele * (1 - theta)):len(R)]
  405. self.p = self.p - 1
  406. # global error
  407. R_g = (C.T.dot(q) + G.dot(u) + H.dot(U_hat[1:-1])).dot(1 - z_hat)
  408. return np.abs(np.sum(R) + np.sum(R_g)), refine_index + 1
  409. def adaptive(self):
  410. x = np.linspace(0, 1, self.n_ele + 1)
  411. index = []
  412. U, hat_U = self.solve_local(index, x)
  413. U_adjoint, hat_adjoint = self.solve_adjoint(index, x, U, hat_U)
  414. X = np.zeros(self.p * self.n_ele)
  415. dx = np.zeros(self.n_ele + len(index))
  416. for i in range(1, self.n_ele + 1):
  417. x_i, dx_i, x_n = self.mesh(i, index, x)
  418. X[(i - 1) * self.p:(i - 1) * self.p + self.p] = x_i
  419. dx[i - 1] = dx_i
  420. numAdaptive = 28
  421. trueError = np.zeros((2, numAdaptive))
  422. estError = np.zeros((2, numAdaptive))
  423. for i in np.arange(numAdaptive):
  424. est_error, index = self.residual(
  425. U, hat_U, U_adjoint, hat_adjoint, dx, index, x)
  426. index = index.tolist()
  427. U, hat_U = self.solve_local(index, x)
  428. U_adjoint, hat_adjoint = self.solve_adjoint(index, x, U, hat_U)
  429. self.p = self.p - 1
  430. X = np.zeros(self.p * (self.n_ele + len(index)))
  431. dx = np.zeros(self.n_ele + len(index))
  432. for j in range(1, self.n_ele + len(index) + 1):
  433. x_i, dx_i, x_n = self.mesh(j, index, x)
  434. X[(j - 1) * self.p:(j - 1) * self.p + self.p] = x_i
  435. dx[j - 1] = dx_i
  436. x = x_n
  437. estError[0, i] = self.n_ele
  438. estError[1, i] = est_error
  439. self.n_ele = self.n_ele + len(index)
  440. # U_1d = np.zeros(len(U))
  441. # for j in np.arange(len(U)):
  442. # U_1d[j] = U[j]
  443. # Unum = np.array([])
  444. # Xnum = np.array([])
  445. # Qnum = np.array([])
  446. # for j in range(1, self.n_ele + 1):
  447. # # Gauss quadrature
  448. # gorder = 10 * self.p
  449. # xi, wi = np.polynomial.legendre.leggauss(gorder)
  450. # shp, shpx = self.shape(xi, self.p)
  451. # Xnum = np.hstack((Xnum, (X[(j - 1) * self.p + self.p - 1] + X[(j - 1) * self.p]) / 2 + (
  452. # X[(j - 1) * self.p + self.p - 1] - X[(j - 1) * self.p]) / 2 * xi))
  453. # Unum = np.hstack(
  454. # (Unum, shp.T.dot(U_1d[int(len(U) / 2) + (j - 1) * self.p:int(len(U) / 2) + j * self.p])))
  455. # Qnum = np.hstack(
  456. # (Qnum, shp.T.dot(U_1d[int((j - 1) * self.p):j * self.p])))
  457. # if i in [0, 4, 9, 19]:
  458. # plt.plot(Xnum, Unum, '-', color='C3')
  459. # plt.plot(X, U[int(len(U) / 2):len(U)], 'C3.')
  460. # plt.xlabel('$x$', fontsize=17)
  461. # plt.ylabel('$u$', fontsize=17)
  462. # plt.axis([-0.05, 1.05, 0, 1.3])
  463. # plt.grid()
  464. # # plt.savefig('u_test_{}.pdf'.format(i+1))
  465. # plt.show()
  466. # plt.clf()
  467. trueError[0, i] = self.n_ele
  468. trueError[1, i] = np.abs(
  469. U[self.n_ele * self.p - 1] - np.sqrt(self.kappa))
  470. self.p = self.p + 1
  471. self.trueError = trueError
  472. self.estError = estError