Block-Structured AMR Software Framework
AMReX_MLEBABecLap_2D_K.H
Go to the documentation of this file.
1 #ifndef AMREX_MLEBABECLAP_2D_K_H_
2 #define AMREX_MLEBABECLAP_2D_K_H_
3 #include <AMReX_Config.H>
4 #include <AMReX_REAL.H>
5 
7 
8 namespace amrex {
9 
11 void mlebabeclap_adotx_centroid (Box const& box, Array4<Real> const& y,
12  Array4<Real const> const& x, Array4<Real const> const& a,
13  Array4<Real const> const& bX, Array4<Real const> const& bY,
14  Array4<EBCellFlag const> const& flag,
15  Array4<Real const> const& vfrc,
16  Array4<Real const> const& apx, Array4<Real const> const& apy,
17  Array4<Real const> const& fcx, Array4<Real const> const& fcy,
18  Array4<Real const> const& ccent, Array4<Real const> const& ba,
19  Array4<Real const> const& bcent, Array4<Real const> const& beb,
20  Array4<Real const> const& phieb,
21  const int& domlo_x, const int& domlo_y,
22  const int& domhi_x, const int& domhi_y,
23  const bool& on_x_face, const bool& on_y_face,
24  bool is_eb_dirichlet, bool is_eb_inhomog,
25  GpuArray<Real,AMREX_SPACEDIM> const& dxinv,
26  Real alpha, Real beta, int ncomp) noexcept
27 {
28  Real dhx = beta*dxinv[0]*dxinv[0];
29  Real dhy = beta*dxinv[1]*dxinv[1];
30  Real dh = beta*dxinv[0]*dxinv[1];
31 
32  amrex::Loop(box, ncomp, [=] (int i, int j, int k, int n) noexcept
33  {
34  if (flag(i,j,k).isCovered())
35  {
36  y(i,j,k,n) = Real(0.0);
37  }
38  else if (flag(i,j,k).isRegular() &&
39  ((flag(i-1,j ,k).isRegular() && flag(i+1,j ,k).isRegular() &&
40  flag(i ,j-1,k).isRegular() && flag(i ,j+1,k).isRegular()) ))
41  {
42  y(i,j,k,n) = alpha*a(i,j,k)*x(i,j,k,n)
43  - dhx * (bX(i+1,j,k,n)*(x(i+1,j,k,n) - x(i ,j,k,n))
44  - bX(i ,j,k,n)*(x(i ,j,k,n) - x(i-1,j,k,n)))
45  - dhy * (bY(i,j+1,k,n)*(x(i,j+1,k,n) - x(i,j ,k,n))
46  - bY(i,j ,k,n)*(x(i,j ,k,n) - x(i,j-1,k,n)));
47 
48 
49  }
50  else
51  {
52  Real kappa = vfrc(i,j,k);
53  Real apxm = apx(i,j,k);
54  Real apxp = apx(i+1,j,k);
55  Real apym = apy(i,j,k);
56  Real apyp = apy(i,j+1,k);
57 
58  // First get EB-aware slope that doesn't know about extdir
59  bool needs_bdry_stencil = (i <= domlo_x) || (i >= domhi_x) ||
60  (j <= domlo_y) || (j >= domhi_y);
61 
62  // if phi_on_centroid -- A second order least squares fit is used
63  // to approximate the slope on the high and low faces. Note that if
64  // any of the three cells --e.g., (i-1,j), (i,j), or (i-1,j)-- are
65  // cut, then the least squares fit is needed. This is a bit more than
66  // is actually needed for most cases but it will return the correct
67  // value in all cases.
68 
69  Real fxm = bX(i,j,k,n) * (x(i,j,k,n)-x(i-1,j,k,n));
70  if ( (apxm != Real(0.0)) && (vfrc(i,j,k) != Real(1.0) || vfrc(i-1,j,k) != Real(1.0) || vfrc(i+1,j,k) != Real(1.0)) )
71  {
72  Real yloc_on_xface = fcx(i,j,k);
73 
74  if(needs_bdry_stencil) {
75 
76  fxm = grad_x_of_phi_on_centroids_extdir(i,j,k,n,x,phieb,flag,ccent,bcent,vfrc,
77  yloc_on_xface,is_eb_dirichlet,is_eb_inhomog,
78  on_x_face, domlo_x, domhi_x,
79  on_y_face, domlo_y, domhi_y);
80 
81  } else {
82  fxm = grad_x_of_phi_on_centroids(i,j,k,n,x,phieb,flag,ccent,bcent,
83  yloc_on_xface,is_eb_dirichlet,is_eb_inhomog);
84  }
85  fxm *= bX(i,j,k,n);
86  }
87 
88  Real fxp = bX(i+1,j,k,n)*(x(i+1,j,k,n)-x(i,j,k,n));
89  if ( (apxp != Real(0.0)) && (vfrc(i,j,k) != Real(1.0) || vfrc(i+1,j,k) != Real(1.0) || vfrc(i-1,j,k) != Real(1.0)) ) {
90  Real yloc_on_xface = fcx(i+1,j,k,0);
91  if(needs_bdry_stencil) {
92  fxp = grad_x_of_phi_on_centroids_extdir(i+1,j,k,n,x,phieb,flag,ccent,bcent,vfrc,
93  yloc_on_xface,is_eb_dirichlet,is_eb_inhomog,
94  on_x_face, domlo_x, domhi_x,
95  on_y_face, domlo_y, domhi_y);
96 
97  } else {
98  fxp = grad_x_of_phi_on_centroids(i+1,j,k,n,x,phieb,flag,ccent,bcent,
99  yloc_on_xface,is_eb_dirichlet,is_eb_inhomog);
100  }
101  fxp *= bX(i+1,j,k,n);
102 
103  }
104 
105  Real fym = bY(i,j,k,n)*(x(i,j,k,n)-x(i,j-1,k,n));
106  if ( (apym != Real(0.0)) && (vfrc(i,j,k) != Real(1.0) || vfrc(i,j-1,k) != Real(1.0) || vfrc(i,j+1,k) != Real(1.0)) ) {
107  Real xloc_on_yface = fcy(i,j,k,0);
108 
109  if(needs_bdry_stencil) {
110 
111  fym = grad_y_of_phi_on_centroids_extdir(i,j,k,n,x,phieb,flag,ccent,bcent,vfrc,
112  xloc_on_yface,is_eb_dirichlet,is_eb_inhomog,
113  on_x_face, domlo_x, domhi_x,
114  on_y_face, domlo_y, domhi_y);
115 
116  } else {
117  fym = grad_y_of_phi_on_centroids(i,j,k,n,x,phieb,flag,ccent,bcent,
118  xloc_on_yface,is_eb_dirichlet,is_eb_inhomog);
119  }
120  fym *= bY(i,j,k,n);
121  }
122 
123  Real fyp = bY(i,j+1,k,n)*(x(i,j+1,k,n)-x(i,j,k,n));
124  if ( (apyp != Real(0.0)) && (vfrc(i,j,k) != Real(1.0) || vfrc(i,j+1,k) != Real(1.0) || vfrc(i,j-1,k) != Real(1.0)) ) {
125  Real xloc_on_yface = fcy(i,j+1,k,0);
126  if(needs_bdry_stencil) {
127  fyp = grad_y_of_phi_on_centroids_extdir(i,j+1,k,n,x,phieb,flag,ccent,bcent,vfrc,
128  xloc_on_yface,is_eb_dirichlet,is_eb_inhomog,
129  on_x_face, domlo_x, domhi_x,
130  on_y_face, domlo_y, domhi_y);
131 
132  } else {
133  fyp = grad_y_of_phi_on_centroids(i,j+1,k,n,x,phieb,flag,ccent,bcent,
134  xloc_on_yface,is_eb_dirichlet,is_eb_inhomog);
135  }
136  fyp *= bY(i,j+1,k,n);
137  }
138 
139  Real feb = Real(0.0);
140  if (is_eb_dirichlet && flag(i,j,k).isSingleValued())
141  {
142  Real dapx = (apxm-apxp)/dxinv[1];
143  Real dapy = (apym-apyp)/dxinv[0];
144  Real anorm = std::hypot(dapx,dapy);
145  Real anorminv = Real(1.0)/anorm;
146  Real anrmx = dapx * anorminv;
147  Real anrmy = dapy * anorminv;
148 
149 
150  feb = grad_eb_of_phi_on_centroids_extdir(i,j,k,n,x,phieb,flag,ccent,bcent,vfrc,
151  anrmx,anrmy,is_eb_inhomog,
152  on_x_face, domlo_x, domhi_x,
153  on_y_face, domlo_y, domhi_y);
154 
155  feb *= ba(i,j,k) * beb(i,j,k,n);
156  }
157 
158 
159  y(i,j,k,n) = alpha*a(i,j,k)*x(i,j,k,n) + (Real(1.0)/kappa) *
160  (dhx*(apxm*fxm-apxp*fxp) + dhy*(apym*fym-apyp*fyp) - dh*feb);
161  }
162  });
163 }
164 
166 void mlebabeclap_adotx (Box const& box, Array4<Real> const& y,
167  Array4<Real const> const& x, Array4<Real const> const& a,
168  Array4<Real const> const& bX, Array4<Real const> const& bY,
169  Array4<const int> const& ccm, Array4<EBCellFlag const> const& flag,
170  Array4<Real const> const& vfrc, Array4<Real const> const& apx,
171  Array4<Real const> const& apy, Array4<Real const> const& fcx,
172  Array4<Real const> const& fcy, Array4<Real const> const& ba,
173  Array4<Real const> const& bc, Array4<Real const> const& beb,
174  bool is_dirichlet, Array4<Real const> const& phieb,
175  bool is_inhomog, GpuArray<Real,AMREX_SPACEDIM> const& dxinv,
176  Real alpha, Real beta, int ncomp,
177  bool beta_on_centroid, bool phi_on_centroid) noexcept
178 {
179  Real dhx = beta*dxinv[0]*dxinv[0];
180  Real dhy = beta*dxinv[1]*dxinv[1];
181  Real dh = beta*dxinv[0]*dxinv[1];
182 
183 
184  bool beta_on_center = !(beta_on_centroid);
185  bool phi_on_center = !( phi_on_centroid);
186 
187  amrex::Loop(box, ncomp, [=] (int i, int j, int k, int n) noexcept
188  {
189  if (flag(i,j,k).isCovered())
190  {
191  y(i,j,k,n) = Real(0.0);
192  }
193  else if (flag(i,j,k).isRegular())
194  {
195  y(i,j,k,n) = alpha*a(i,j,k)*x(i,j,k,n)
196  - dhx * (bX(i+1,j,k,n)*(x(i+1,j,k,n) - x(i ,j,k,n))
197  - bX(i ,j,k,n)*(x(i ,j,k,n) - x(i-1,j,k,n)))
198  - dhy * (bY(i,j+1,k,n)*(x(i,j+1,k,n) - x(i,j ,k,n))
199  - bY(i,j ,k,n)*(x(i,j ,k,n) - x(i,j-1,k,n)));
200  }
201  else
202  {
203  Real kappa = vfrc(i,j,k);
204  Real apxm = apx(i,j,k);
205  Real apxp = apx(i+1,j,k);
206  Real apym = apy(i,j,k);
207  Real apyp = apy(i,j+1,k);
208 
209  Real fxm = bX(i,j,k,n) * (x(i,j,k,n)-x(i-1,j,k,n));
210  if (apxm != Real(0.0) && apxm != Real(1.0)) {
211  int jj = j + static_cast<int>(std::copysign(Real(1.0),fcx(i,j,k)));
212  Real fracy = (ccm(i-1,jj,k) || ccm(i,jj,k)) ? std::abs(fcx(i,j,k)) : Real(0.0);
213  if (beta_on_center && phi_on_center) {
214  fxm = (Real(1.0)-fracy)*fxm + fracy*bX(i,jj,k,n)*(x(i,jj,k,n)-x(i-1,jj,k,n));
215  } else if (beta_on_centroid && phi_on_center) {
216  fxm = bX(i,j,k,n) * ( (Real(1.0)-fracy)*(x(i, j,k,n)-x(i-1, j,k,n))
217  + fracy *(x(i,jj,k,n)-x(i-1,jj,k,n)) );
218  }
219  }
220 
221  Real fxp = bX(i+1,j,k,n)*(x(i+1,j,k,n)-x(i,j,k,n));
222  if (apxp != Real(0.0) && apxp != Real(1.0)) {
223  int jj = j + static_cast<int>(std::copysign(Real(1.0),fcx(i+1,j,k)));
224  Real fracy = (ccm(i,jj,k) || ccm(i+1,jj,k)) ? std::abs(fcx(i+1,j,k)) : Real(0.0);
225  if (beta_on_center && phi_on_center) {
226  fxp = (Real(1.0)-fracy)*fxp + fracy*bX(i+1,jj,k,n)*(x(i+1,jj,k,n)-x(i,jj,k,n));
227  } else if (beta_on_centroid && phi_on_center) {
228  fxp = bX(i+1,j,k,n) * ( (Real(1.0)-fracy)*(x(i+1, j,k,n)-x(i, j,k,n))
229  + fracy *(x(i+1,jj,k,n)-x(i,jj,k,n)) );
230  }
231  }
232 
233  Real fym = bY(i,j,k,n)*(x(i,j,k,n)-x(i,j-1,k,n));
234  if (apym != Real(0.0) && apym != Real(1.0)) {
235  int ii = i + static_cast<int>(std::copysign(Real(1.0),fcy(i,j,k)));
236  Real fracx = (ccm(ii,j-1,k) || ccm(ii,j,k)) ? std::abs(fcy(i,j,k)) : Real(0.0);
237  if (beta_on_center && phi_on_center) {
238  fym = (Real(1.0)-fracx)*fym + fracx*bY(ii,j,k,n)*(x(ii,j,k,n)-x(ii,j-1,k,n));
239  } else if (beta_on_centroid && phi_on_center) {
240  fym = bY(i,j,k,n) * ( (Real(1.0)-fracx)*(x( i,j,k,n)-x( i,j-1,k,n))
241  + fracx *(x(ii,j,k,n)-x(ii,j-1,k,n)) );
242  }
243  }
244 
245  Real fyp = bY(i,j+1,k,n)*(x(i,j+1,k,n)-x(i,j,k,n));
246  if (apyp != Real(0.0) && apyp != Real(1.0)) {
247  int ii = i + static_cast<int>(std::copysign(Real(1.0),fcy(i,j+1,k)));
248  Real fracx = (ccm(ii,j,k) || ccm(ii,j+1,k)) ? std::abs(fcy(i,j+1,k)) : Real(0.0);
249  if (beta_on_center && phi_on_center) {
250  fyp = (Real(1.0)-fracx)*fyp + fracx*bY(ii,j+1,k,n)*(x(ii,j+1,k,n)-x(ii,j,k,n));
251  } else if (beta_on_centroid && phi_on_center) {
252  fyp = bY(i,j+1,k,n) * ( (Real(1.0)-fracx)*(x( i,j+1,k,n)-x( i,j,k,n))
253  + fracx *(x(ii,j+1,k,n)-x(ii,j,k,n)) );
254  }
255  }
256 
257  Real feb = Real(0.0);
258  if (is_dirichlet) {
259  Real dapx = (apxm-apxp)/dxinv[1];
260  Real dapy = (apym-apyp)/dxinv[0];
261  Real anorm = std::hypot(dapx,dapy);
262  Real anorminv = Real(1.0)/anorm;
263  Real anrmx = dapx * anorminv;
264  Real anrmy = dapy * anorminv;
265 
266  Real phib = is_inhomog ? phieb(i,j,k,n) : Real(0.0);
267 
268  Real bctx = bc(i,j,k,0);
269  Real bcty = bc(i,j,k,1);
270  Real dx_eb = get_dx_eb(kappa);
271 
272  Real dg, gx, gy, sx, sy;
273  if (std::abs(anrmx) > std::abs(anrmy)) {
274  dg = dx_eb / std::abs(anrmx);
275  } else {
276  dg = dx_eb / std::abs(anrmy);
277  }
278  gx = (bctx - dg*anrmx);
279  gy = (bcty - dg*anrmy);
280  sx = std::copysign(Real(1.0),anrmx);
281  sy = std::copysign(Real(1.0),anrmy);
282 
283  int ii = i - static_cast<int>(sx);
284  int jj = j - static_cast<int>(sy);
285 
286  Real phig = (Real(1.0) + gx*sx + gy*sy + gx*gy*sx*sy) * x(i ,j ,k,n)
287  + ( - gx*sx - gx*gy*sx*sy) * x(ii,j ,k,n)
288  + ( - gy*sy - gx*gy*sx*sy) * x(i ,jj,k,n)
289  + ( + gx*gy*sx*sy) * x(ii,jj,k,n) ;
290 
291  Real dphidn = (phib-phig) / dg;
292 
293  feb = dphidn * ba(i,j,k) * beb(i,j,k,n);
294  }
295 
296  y(i,j,k,n) = alpha*a(i,j,k)*x(i,j,k,n) + (Real(1.0)/kappa) *
297  (dhx*(apxm*fxm-apxp*fxp) +
298  dhy*(apym*fym-apyp*fyp) - dh*feb);
299  }
300  });
301 }
302 
304 void mlebabeclap_ebflux (int i, int j, int k, int n,
305  Array4<Real> const& feb,
306  Array4<Real const> const& x,
307  Array4<EBCellFlag const> const& flag,
308  Array4<Real const> const& vfrc,
309  Array4<Real const> const& apx,
310  Array4<Real const> const& apy,
311  Array4<Real const> const& bc,
312  Array4<Real const> const& beb,
313  Array4<Real const> const& phieb,
314  bool is_inhomog,
315  GpuArray<Real,AMREX_SPACEDIM> const& dxinv) noexcept
316 {
317 
318  if (!flag(i,j,k).isSingleValued())
319  {
320  feb(i,j,k,n) = Real(0.0);
321  }
322  else
323  {
324  Real kappa = vfrc(i,j,k);
325  Real apxm = apx(i,j,k);
326  Real apxp = apx(i+1,j,k);
327  Real apym = apy(i,j,k);
328  Real apyp = apy(i,j+1,k);
329 
330  Real dapx = (apxm-apxp)/dxinv[1];
331  Real dapy = (apym-apyp)/dxinv[0];
332  Real anorm = std::hypot(dapx,dapy);
333  Real anorminv = Real(1.0)/anorm;
334  Real anrmx = dapx * anorminv;
335  Real anrmy = dapy * anorminv;
336  const Real bareascaling = std::sqrt( (anrmx/dxinv[0])*(anrmx/dxinv[0]) +
337  (anrmy/dxinv[1])*(anrmy/dxinv[1]) );
338 
339 
340  Real phib = is_inhomog ? phieb(i,j,k,n) : Real(0.0);
341 
342  Real bctx = bc(i,j,k,0);
343  Real bcty = bc(i,j,k,1);
344  Real dx_eb = get_dx_eb(kappa);
345 
346  Real dg, gx, gy, sx, sy;
347  if (std::abs(anrmx) > std::abs(anrmy)) {
348  dg = dx_eb / std::abs(anrmx);
349  } else {
350  dg = dx_eb / std::abs(anrmy);
351  }
352  gx = bctx - dg*anrmx;
353  gy = bcty - dg*anrmy;
354  sx = std::copysign(Real(1.0),anrmx);
355  sy = std::copysign(Real(1.0),anrmy);
356 
357  int ii = i - static_cast<int>(sx);
358  int jj = j - static_cast<int>(sy);
359 
360  Real phig = (Real(1.0) + gx*sx + gy*sy + gx*gy*sx*sy) * x(i ,j ,k,n)
361  + ( - gx*sx - gx*gy*sx*sy) * x(ii,j ,k,n)
362  + ( - gy*sy - gx*gy*sx*sy) * x(i ,jj,k,n)
363  + ( + gx*gy*sx*sy) * x(ii,jj,k,n) ;
364 
365  Real dphidn = (phib-phig)/(dg * bareascaling);
366  feb(i,j,k,n) = -beb(i,j,k,n) * dphidn;
367  }
368 }
369 
371 void mlebabeclap_gsrb (Box const& box,
372  Array4<Real> const& phi, Array4<Real const> const& rhs,
373  Real alpha, Array4<Real const> const& a,
374  Real dhx, Real dhy, Real dh,
376  Array4<Real const> const& bX, Array4<Real const> const& bY,
377  Array4<int const> const& m0, Array4<int const> const& m2,
378  Array4<int const> const& m1, Array4<int const> const& m3,
379  Array4<Real const> const& f0, Array4<Real const> const& f2,
380  Array4<Real const> const& f1, Array4<Real const> const& f3,
381  Array4<const int> const& ccm, Array4<Real const> const& beb,
382  EBData const& ebdata,
383  bool is_dirichlet, bool beta_on_centroid, bool phi_on_centroid,
384  Box const& vbox, int redblack, int ncomp) noexcept
385 {
386  const auto vlo = amrex::lbound(vbox);
387  const auto vhi = amrex::ubound(vbox);
388 
389  amrex::Loop(box, ncomp, [=] (int i, int j, int k, int n) noexcept
390  {
391  if ((i+j+k+redblack) % 2 == 0)
392  {
393  auto const flag = ebdata.get<EBData_t::cellflag>(i,j,k);
394  if (flag.isCovered())
395  {
396  phi(i,j,k,n) = Real(0.0);
397  }
398  else
399  {
400  Real cf0 = (i == vlo.x && m0(vlo.x-1,j,k) > 0)
401  ? f0(vlo.x,j,k,n) : Real(0.0);
402  Real cf1 = (j == vlo.y && m1(i,vlo.y-1,k) > 0)
403  ? f1(i,vlo.y,k,n) : Real(0.0);
404  Real cf2 = (i == vhi.x && m2(vhi.x+1,j,k) > 0)
405  ? f2(vhi.x,j,k,n) : Real(0.0);
406  Real cf3 = (j == vhi.y && m3(i,vhi.y+1,k) > 0)
407  ? f3(i,vhi.y,k,n) : Real(0.0);
408 
409  if (flag.isRegular())
410  {
411  Real gamma = alpha*a(i,j,k)
412  + dhx * (bX(i+1,j,k,n) + bX(i,j,k,n))
413  + dhy * (bY(i,j+1,k,n) + bY(i,j,k,n));
414 
415  Real rho = dhx * (bX(i+1,j,k,n)*phi(i+1,j,k,n)
416  + bX(i ,j,k,n)*phi(i-1,j,k,n))
417  + dhy * (bY(i,j+1,k,n)*phi(i,j+1,k,n)
418  + bY(i,j ,k,n)*phi(i,j-1,k,n));
419 
420  Real delta = dhx*(bX(i,j,k,n)*cf0 + bX(i+1,j,k,n)*cf2)
421  + dhy*(bY(i,j,k,n)*cf1 + bY(i,j+1,k,n)*cf3);
422 
423  Real res = rhs(i,j,k,n) - (gamma*phi(i,j,k,n) - rho);
424  phi(i,j,k,n) += res/(gamma-delta);
425  }
426  else
427  {
428  Real kappa = ebdata.get<EBData_t::volfrac>(i,j,k);
429  Real apxm = ebdata.get<EBData_t::apx>(i ,j ,k);
430  Real apxp = ebdata.get<EBData_t::apx>(i+1,j ,k);
431  Real apym = ebdata.get<EBData_t::apy>(i ,j ,k);
432  Real apyp = ebdata.get<EBData_t::apy>(i ,j+1,k);
433 
434  Real fxm = -bX(i,j,k,n)*phi(i-1,j,k,n);
435  Real oxm = -bX(i,j,k,n)*cf0;
436  Real sxm = bX(i,j,k,n);
437  if (apxm != Real(0.0) && apxm != Real(1.0)) {
438  Real fcx = ebdata.get<EBData_t::fcx>(i,j,k);
439  int jj = j + static_cast<int>(std::copysign(Real(1.0),fcx));
440  Real fracy = (ccm(i-1,jj,k) || ccm(i,jj,k)) ? std::abs(fcx) : Real(0.0);
441  if (!beta_on_centroid && !phi_on_centroid)
442  {
443  fxm = (Real(1.0)-fracy)*fxm +
444  fracy *bX(i,jj,k,n)*(phi(i,jj,k,n)-phi(i-1,jj,k,n));
445  }
446  else if (beta_on_centroid && !phi_on_centroid)
447  {
448  fxm = (Real(1.0)-fracy)*( -phi(i-1,j,k,n)) +
449  fracy *(phi(i,jj,k,n)-phi(i-1,jj,k,n));
450  fxm *= bX(i,j,k,n);
451  }
452  oxm = Real(0.0);
453  sxm = (Real(1.0)-fracy)*sxm;
454  }
455 
456  Real fxp = bX(i+1,j,k,n)*phi(i+1,j,k,n);
457  Real oxp = bX(i+1,j,k,n)*cf2;
458  Real sxp = -bX(i+1,j,k,n);
459  if (apxp != Real(0.0) && apxp != Real(1.0)) {
460  Real fcx = ebdata.get<EBData_t::fcx>(i+1,j,k);
461  int jj = j + static_cast<int>(std::copysign(Real(1.0),fcx));
462  Real fracy = (ccm(i,jj,k) || ccm(i+1,jj,k)) ? std::abs(fcx) : Real(0.0);
463  if (!beta_on_centroid && !phi_on_centroid)
464  {
465  fxp = (Real(1.0)-fracy)*fxp +
466  fracy *bX(i+1,jj,k,n)*(phi(i+1,jj,k,n)-phi(i,jj,k,n));
467  }
468  else if (beta_on_centroid && !phi_on_centroid)
469  {
470  fxp = (Real(1.0)-fracy)*(phi(i+1,j,k,n) ) +
471  fracy *(phi(i+1,jj,k,n)-phi(i,jj,k,n));
472  fxp *= bX(i+1,j,k,n);
473  }
474  oxp = Real(0.0);
475  sxp = (Real(1.0)-fracy)*sxp;
476  }
477 
478  Real fym = -bY(i,j,k,n)*phi(i,j-1,k,n);
479  Real oym = -bY(i,j,k,n)*cf1;
480  Real sym = bY(i,j,k,n);
481  if (apym != Real(0.0) && apym != Real(1.0)) {
482  Real fcy = ebdata.get<EBData_t::fcy>(i,j,k);
483  int ii = i + static_cast<int>(std::copysign(Real(1.0),fcy));
484  Real fracx = (ccm(ii,j-1,k) || ccm(ii,j,k)) ? std::abs(fcy) : Real(0.0);
485  if (!beta_on_centroid && !phi_on_centroid)
486  {
487  fym = (Real(1.0)-fracx)*fym +
488  fracx *bY(ii,j,k,n)*(phi(ii,j,k,n)-phi(ii,j-1,k,n));
489  }
490  else if (beta_on_centroid && !phi_on_centroid)
491  {
492  fym = (Real(1.0)-fracx)*( -phi( i,j-1,k,n)) +
493  fracx *(phi(ii,j,k,n)-phi(ii,j-1,k,n));
494  fym *= bY(i,j,k,n);
495  }
496 
497  oym = Real(0.0);
498  sym = (Real(1.0)-fracx)*sym;
499  }
500 
501  Real fyp = bY(i,j+1,k,n)*phi(i,j+1,k,n);
502  Real oyp = bY(i,j+1,k,n)*cf3;
503  Real syp = -bY(i,j+1,k,n);
504  if (apyp != Real(0.0) && apyp != Real(1.0)) {
505  Real fcy = ebdata.get<EBData_t::fcy>(i,j+1,k);
506  int ii = i + static_cast<int>(std::copysign(Real(1.0),fcy));
507  Real fracx = (ccm(ii,j,k) || ccm(ii,j+1,k)) ? std::abs(fcy) : Real(0.0);
508  if (!beta_on_centroid && !phi_on_centroid)
509  {
510  fyp = (Real(1.0)-fracx)*fyp +
511  fracx*bY(ii,j+1,k,n)*(phi(ii,j+1,k,n)-phi(ii,j,k,n));
512  }
513  else if (beta_on_centroid && !phi_on_centroid)
514  {
515  fyp = (Real(1.0)-fracx)*(phi(i,j+1,k,n) )+
516  fracx *(phi(ii,j+1,k,n)-phi(ii,j,k,n));
517  fyp *= bY(i,j+1,k,n);
518  }
519  oyp = Real(0.0);
520  syp = (Real(1.0)-fracx)*syp;
521  }
522 
523  Real vfrcinv = (Real(1.0)/kappa);
524  Real gamma = alpha*a(i,j,k) + vfrcinv *
525  (dhx*(apxm*sxm-apxp*sxp) +
526  dhy*(apym*sym-apyp*syp));
527  Real rho = -vfrcinv *
528  (dhx*(apxm*fxm-apxp*fxp) +
529  dhy*(apym*fym-apyp*fyp));
530 
531  Real delta = -vfrcinv *
532  (dhx*(apxm*oxm-apxp*oxp) +
533  dhy*(apym*oym-apyp*oyp));
534 
535  if (is_dirichlet) {
536  Real dapx = (apxm-apxp)*dx[1];
537  Real dapy = (apym-apyp)*dx[0];
538  Real anorm = std::hypot(dapx,dapy);
539  Real anorminv = Real(1.0)/anorm;
540  Real anrmx = dapx * anorminv;
541  Real anrmy = dapy * anorminv;
542 
543  Real bctx = ebdata.get<EBData_t::bndrycent>(i,j,k,0);
544  Real bcty = ebdata.get<EBData_t::bndrycent>(i,j,k,1);
545  Real dx_eb = get_dx_eb(kappa);
546 
547  Real dg, gx, gy, sx, sy;
548  if (std::abs(anrmx) > std::abs(anrmy)) {
549  dg = dx_eb / std::abs(anrmx);
550  } else {
551  dg = dx_eb / std::abs(anrmy);
552  }
553  gx = bctx - dg*anrmx;
554  gy = bcty - dg*anrmy;
555  sx = std::copysign(Real(1.0),anrmx);
556  sy = std::copysign(Real(1.0),anrmy);
557 
558  int ii = i - static_cast<int>(sx);
559  int jj = j - static_cast<int>(sy);
560 
561  Real phig_gamma = (Real(1.0) + gx*sx + gy*sy + gx*gy*sx*sy);
562  Real phig = ( - gx*sx - gx*gy*sx*sy) * phi(ii,j ,k,n)
563  + ( - gy*sy - gx*gy*sx*sy) * phi(i ,jj,k,n)
564  + ( + gx*gy*sx*sy) * phi(ii,jj,k,n);
565 
566  // In gsrb we are always in residual-correction form so phib = 0
567  Real dphidn = ( -phig)/dg;
568 
569  Real ba = ebdata.get<EBData_t::bndryarea>(i,j,k);
570 
571  Real feb = dphidn * ba * beb(i,j,k,n);
572  rho += -vfrcinv*(-dh)*feb;
573 
574  Real feb_gamma = -phig_gamma/dg * ba * beb(i,j,k,n);
575  gamma += vfrcinv*(-dh)*feb_gamma;
576  }
577 
578  Real res = rhs(i,j,k,n) - (gamma*phi(i,j,k,n) - rho);
579  phi(i,j,k,n) += res/(gamma-delta);
580  }
581  }
582  }
583  });
584 }
585 
587 void mlebabeclap_flux_x (Box const& box, Array4<Real> const& fx, Array4<Real const> const& apx,
588  Array4<Real const> const& fcx, Array4<Real const> const& sol,
589  Array4<Real const> const& bX, Array4<int const> const& ccm,
590  Real dhx, int face_only, int ncomp, Box const& xbox,
591  bool beta_on_centroid, bool phi_on_centroid) noexcept
592 {
593  int lof = xbox.smallEnd(0);
594  int hif = xbox.bigEnd(0);
595  amrex::LoopConcurrent(box, ncomp, [=] (int i, int j, int k, int n) noexcept
596  {
597  if (!face_only || lof == i || hif == i) {
598  if (apx(i,j,k) == Real(0.0)) {
599  fx(i,j,k,n) = Real(0.0);
600  } else if (apx(i,j,k) == Real(1.0)) {
601  fx(i,j,k,n) = -dhx*bX(i,j,k,n)*(sol(i,j,k,n)-sol(i-1,j,k,n));
602  } else {
603  Real fxm = bX(i,j,k,n)*(sol(i,j,k,n)-sol(i-1,j,k,n));
604  int jj = j + static_cast<int>(std::copysign(Real(1.0),fcx(i,j,k)));
605  Real fracy = (ccm(i-1,jj,k) || ccm(i,jj,k)) ? std::abs(fcx(i,j,k)) : Real(0.0);
606  if (!beta_on_centroid && !phi_on_centroid)
607  {
608  fxm = (Real(1.0)-fracy)*fxm + fracy*bX(i,jj,k,n)*(sol(i,jj,k,n)-sol(i-1,jj,k,n));
609  } else if (beta_on_centroid && !phi_on_centroid) {
610  fxm = bX(i,j,k,n) * ( (Real(1.0)-fracy)*(sol(i, j,k,n)-sol(i-1, j,k,n)) +
611  fracy *(sol(i,jj,k,n)-sol(i-1,jj,k,n)) );
612  }
613  fx(i,j,k,n) = -fxm*dhx;
614  }
615  }
616  });
617 }
618 
620 void mlebabeclap_flux_y (Box const& box, Array4<Real> const& fy, Array4<Real const> const& apy,
621  Array4<Real const> const& fcy, Array4<Real const> const& sol,
622  Array4<Real const> const& bY, Array4<int const> const& ccm,
623  Real dhy, int face_only, int ncomp, Box const& ybox,
624  bool beta_on_centroid, bool phi_on_centroid) noexcept
625 {
626  int lof = ybox.smallEnd(1);
627  int hif = ybox.bigEnd(1);
628  amrex::LoopConcurrent(box, ncomp, [=] (int i, int j, int k, int n) noexcept
629  {
630  if (!face_only || lof == j || hif == j) {
631  if (apy(i,j,k) == Real(0.0)) {
632  fy(i,j,k,n) = Real(0.0);
633  } else if (apy(i,j,k) == Real(1.0)) {
634  fy(i,j,k,n) = -dhy*bY(i,j,k,n)*(sol(i,j,k,n)-sol(i,j-1,k,n));
635  } else {
636  Real fym = bY(i,j,k,n)*(sol(i,j,k,n)-sol(i,j-1,k,n));
637  int ii = i + static_cast<int>(std::copysign(Real(1.0),fcy(i,j,k)));
638  Real fracx = (ccm(ii,j-1,k) || ccm(ii,j,k)) ? std::abs(fcy(i,j,k)) : Real(0.0);
639  if (!beta_on_centroid && !phi_on_centroid)
640  {
641  fym = (Real(1.0)-fracx)*fym + fracx*bY(ii,j,k,n)*(sol(ii,j,k,n)-sol(ii,j-1,k,n));
642  } else if (beta_on_centroid && !phi_on_centroid) {
643  fym = bY(i,j,k,n) * ( (Real(1.0)-fracx)*(sol( i,j,k,n)-sol( i,j-1,k,n)) +
644  fracx *(sol(ii,j,k,n)-sol(ii,j-1,k,n)) );
645  }
646  fy(i,j,k,n) = -fym*dhy;
647  }
648  }
649  });
650 }
651 
653 void mlebabeclap_flux_x_0 (Box const& box, Array4<Real> const& fx, Array4<Real const> const& apx,
654  Array4<Real const> const& sol, Array4<Real const> const& bX,
655  Real dhx, int face_only, int ncomp, Box const& xbox) noexcept
656 {
657  int lof = xbox.smallEnd(0);
658  int hif = xbox.bigEnd(0);
659  amrex::LoopConcurrent(box, ncomp, [=] (int i, int j, int k, int n) noexcept
660  {
661  if (!face_only || lof == i || hif == i) {
662  if (apx(i,j,k) == Real(0.0)) {
663  fx(i,j,k,n) = Real(0.0);
664  } else {
665  fx(i,j,k,n) = -dhx*bX(i,j,k,n)*(sol(i,j,k,n)-sol(i-1,j,k,n));
666  }
667  }
668  });
669 }
670 
672 void mlebabeclap_flux_y_0 (Box const& box, Array4<Real> const& fy, Array4<Real const> const& apy,
673  Array4<Real const> const& sol, Array4<Real const> const& bY,
674  Real dhy, int face_only, int ncomp, Box const& ybox) noexcept
675 {
676  int lof = ybox.smallEnd(1);
677  int hif = ybox.bigEnd(1);
678  amrex::LoopConcurrent(box, ncomp, [=] (int i, int j, int k, int n) noexcept
679  {
680  if (!face_only || lof == j || hif == j) {
681  if (apy(i,j,k) == Real(0.0)) {
682  fy(i,j,k,n) = Real(0.0);
683  } else {
684  fy(i,j,k,n) = -dhy*bY(i,j,k,n)*(sol(i,j,k,n)-sol(i,j-1,k,n));
685  }
686  }
687  });
688 }
689 
691 void mlebabeclap_grad_x (Box const& box, Array4<Real> const& gx, Array4<Real const> const& sol,
692  Array4<Real const> const& apx, Array4<Real const> const& fcx,
693  Array4<int const> const& ccm,
694  Real dxi, int ncomp, bool phi_on_centroid) noexcept
695 {
696  amrex::LoopConcurrent(box, ncomp, [=] (int i, int j, int k, int n) noexcept
697  {
698  if (apx(i,j,k) == Real(0.0)) {
699  gx(i,j,k,n) = Real(0.0);
700  } else if (apx(i,j,k) == Real(1.0)) {
701  gx(i,j,k,n) = dxi*(sol(i,j,k,n)-sol(i-1,j,k,n));
702  } else {
703  Real gxm = (sol(i,j,k,n)-sol(i-1,j,k,n));
704  int jj = j + static_cast<int>(std::copysign(Real(1.0),fcx(i,j,k)));
705  Real fracy = (ccm(i-1,jj,k) || ccm(i,jj,k)) ? std::abs(fcx(i,j,k)) : Real(0.0);
706  if (!phi_on_centroid) {
707  gxm = (Real(1.0)-fracy)*gxm + fracy*(sol(i,jj,k,n)-sol(i-1,jj,k,n));
708  }
709  gx(i,j,k,n) = gxm*dxi;
710  }
711  });
712 }
713 
715 void mlebabeclap_grad_y (Box const& box, Array4<Real> const& gy, Array4<Real const> const& sol,
716  Array4<Real const> const& apy, Array4<Real const> const& fcy,
717  Array4<int const> const& ccm,
718  Real dyi, int ncomp, bool phi_on_centroid) noexcept
719 {
720  amrex::LoopConcurrent(box, ncomp, [=] (int i, int j, int k, int n) noexcept
721  {
722  if (apy(i,j,k) == Real(0.0)) {
723  gy(i,j,k,n) = Real(0.0);
724  } else if (apy(i,j,k) == Real(1.0)) {
725  gy(i,j,k,n) = dyi*(sol(i,j,k,n)-sol(i,j-1,k,n));
726  } else {
727  Real gym = (sol(i,j,k,n)-sol(i,j-1,k,n));
728  int ii = i + static_cast<int>(std::copysign(Real(1.0),fcy(i,j,k)));
729  Real fracx = (ccm(ii,j-1,k) || ccm(ii,j,k)) ? std::abs(fcy(i,j,k)) : Real(0.0);
730  if (!phi_on_centroid) {
731  gym = (Real(1.0)-fracx)*gym + fracx*(sol(ii,j,k,n)-sol(ii,j-1,k,n));
732  }
733  gy(i,j,k,n) = gym*dyi;
734  }
735  });
736 }
737 
739 void mlebabeclap_grad_x_0 (Box const& box, Array4<Real> const& gx, Array4<Real const> const& sol,
740  Array4<Real const> const& apx, Real dxi, int ncomp) noexcept
741 {
742  amrex::LoopConcurrent(box, ncomp, [=] (int i, int j, int k, int n) noexcept
743  {
744  if (apx(i,j,k) == Real(0.0)) {
745  gx(i,j,k,n) = Real(0.0);
746  } else {
747  gx(i,j,k,n) = dxi*(sol(i,j,k,n)-sol(i-1,j,k,n));
748  }
749  });
750 }
751 
753 void mlebabeclap_grad_y_0 (Box const& box, Array4<Real> const& gy, Array4<Real const> const& sol,
754  Array4<Real const> const& apy, Real dyi, int ncomp) noexcept
755 {
756  amrex::LoopConcurrent(box, ncomp, [=] (int i, int j, int k, int n) noexcept
757  {
758  if (apy(i,j,k) == Real(0.0)) {
759  gy(i,j,k,n) = Real(0.0);
760  } else {
761  gy(i,j,k,n) = dyi*(sol(i,j,k,n)-sol(i,j-1,k,n));
762  }
763  });
764 }
765 
767 void mlebabeclap_normalize (Box const& box, Array4<Real> const& phi,
768  Real alpha, Array4<Real const> const& a,
769  Real dhx, Real dhy, Real dh,
771  Array4<Real const> const& bX, Array4<Real const> const& bY,
772  Array4<const int> const& ccm, Array4<EBCellFlag const> const& flag,
773  Array4<Real const> const& vfrc,
774  Array4<Real const> const& apx, Array4<Real const> const& apy,
775  Array4<Real const> const& fcx, Array4<Real const> const& fcy,
776  Array4<Real const> const& ba, Array4<Real const> const& bc,
777  Array4<Real const> const& beb,
778  bool is_dirichlet, bool beta_on_centroid, int ncomp) noexcept
779 {
780  amrex::Loop(box, ncomp, [=] (int i, int j, int k, int n) noexcept
781  {
782  if (flag(i,j,k).isRegular())
783  {
784  phi(i,j,k,n) /= alpha*a(i,j,k) + dhx*(bX(i,j,k,n) + bX(i+1,j,k,n))
785  + dhy*(bY(i,j,k,n) + bY(i,j+1,k,n));
786  }
787  else if (flag(i,j,k).isSingleValued())
788  {
789  Real kappa = vfrc(i,j,k);
790  Real apxm = apx(i,j,k);
791  Real apxp = apx(i+1,j,k);
792  Real apym = apy(i,j,k);
793  Real apyp = apy(i,j+1,k);
794 
795  Real sxm = bX(i,j,k,n);
796  if (apxm != Real(0.0) && apxm != Real(1.0) && !beta_on_centroid) {
797  int jj = j + static_cast<int>(std::copysign(Real(1.0),fcx(i,j,k)));
798  Real fracy = (ccm(i-1,jj,k) || ccm(i,jj,k))
799  ? std::abs(fcx(i,j,k)) : Real(0.0);
800  sxm = (Real(1.0)-fracy)*sxm;
801  }
802 
803  Real sxp = -bX(i+1,j,k,n);
804  if (apxp != Real(0.0) && apxp != Real(1.0) && !beta_on_centroid) {
805  int jj = j + static_cast<int>(std::copysign(Real(1.0),fcx(i+1,j,k)));
806  Real fracy = (ccm(i,jj,k) || ccm(i+1,jj,k))
807  ? std::abs(fcx(i+1,j,k)) : Real(0.0);
808  sxp = (Real(1.0)-fracy)*sxp;
809  }
810 
811  Real sym = bY(i,j,k,n);
812  if (apym != Real(0.0) && apym != Real(1.0) && !beta_on_centroid) {
813  int ii = i + static_cast<int>(std::copysign(Real(1.0),fcy(i,j,k)));
814  Real fracx = (ccm(ii,j-1,k) || ccm(ii,j,k))
815  ? std::abs(fcy(i,j,k)) : Real(0.0);
816  sym = (Real(1.0)-fracx)*sym;
817  }
818 
819  Real syp = -bY(i,j+1,k,n);
820  if (apyp != Real(0.0) && apyp != Real(1.0) && !beta_on_centroid) {
821  int ii = i + static_cast<int>(std::copysign(Real(1.0),fcy(i,j+1,k)));
822  Real fracx = (ccm(ii,j,k) || ccm(ii,j+1,k))
823  ? std::abs(fcy(i,j+1,k)) : Real(0.0);
824  syp = (Real(1.0)-fracx)*syp;
825  }
826 
827  Real vfrcinv = (Real(1.0)/kappa);
828  Real gamma = alpha*a(i,j,k) + vfrcinv *
829  (dhx*(apxm*sxm-apxp*sxp) +
830  dhy*(apym*sym-apyp*syp));
831 
832  if (is_dirichlet) {
833  Real dapx = (apxm-apxp)*dx[1];
834  Real dapy = (apym-apyp)*dx[0];
835  Real anorm = std::hypot(dapx,dapy);
836  Real anorminv = Real(1.0)/anorm;
837  Real anrmx = dapx * anorminv;
838  Real anrmy = dapy * anorminv;
839 
840  Real bctx = bc(i,j,k,0);
841  Real bcty = bc(i,j,k,1);
842  Real dx_eb = get_dx_eb(vfrc(i,j,k));
843 
844  Real dg, gx, gy, sx, sy;
845  if (std::abs(anrmx) > std::abs(anrmy)) {
846  dg = dx_eb / std::abs(anrmx);
847  } else {
848  dg = dx_eb / std::abs(anrmy);
849  }
850  gx = bctx - dg*anrmx;
851  gy = bcty - dg*anrmy;
852  sx = std::copysign(Real(1.0),anrmx);
853  sy = std::copysign(Real(1.0),anrmy);
854 
855  Real phig_gamma = (Real(1.0) + gx*sx + gy*sy + gx*gy*sx*sy);
856  Real feb_gamma = -phig_gamma/dg * ba(i,j,k) * beb(i,j,k,n);
857  gamma += vfrcinv*(-dh)*feb_gamma;
858  }
859 
860  phi(i,j,k,n) /= gamma;
861  }
862  });
863 }
864 
865 }
866 #endif
#define AMREX_FORCE_INLINE
Definition: AMReX_Extension.H:119
#define AMREX_GPU_HOST_DEVICE
Definition: AMReX_GpuQualifiers.H:20
Definition: AMReX_Amr.cpp:49
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlebabeclap_grad_x_0(Box const &box, Array4< Real > const &gx, Array4< Real const > const &sol, Array4< Real const > const &apx, Real dxi, int ncomp) noexcept
Definition: AMReX_MLEBABecLap_2D_K.H:739
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real grad_y_of_phi_on_centroids(int i, int j, int k, int n, Array4< Real const > const &phi, Array4< Real const > const &phieb, Array4< EBCellFlag const > const &flag, Array4< Real const > const &ccent, Array4< Real const > const &bcent, Real &xloc_on_yface, bool is_eb_dirichlet, bool is_eb_inhomog)
Definition: AMReX_EB_LeastSquares_2D_K.H:331
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlebabeclap_grad_y_0(Box const &box, Array4< Real > const &gy, Array4< Real const > const &sol, Array4< Real const > const &apy, Real dyi, int ncomp) noexcept
Definition: AMReX_MLEBABecLap_2D_K.H:753
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlebabeclap_adotx_centroid(Box const &box, Array4< Real > const &y, Array4< Real const > const &x, Array4< Real const > const &a, Array4< Real const > const &bX, Array4< Real const > const &bY, Array4< EBCellFlag const > const &flag, Array4< Real const > const &vfrc, Array4< Real const > const &apx, Array4< Real const > const &apy, Array4< Real const > const &fcx, Array4< Real const > const &fcy, Array4< Real const > const &ccent, Array4< Real const > const &ba, Array4< Real const > const &bcent, Array4< Real const > const &beb, Array4< Real const > const &phieb, const int &domlo_x, const int &domlo_y, const int &domhi_x, const int &domhi_y, const bool &on_x_face, const bool &on_y_face, bool is_eb_dirichlet, bool is_eb_inhomog, GpuArray< Real, AMREX_SPACEDIM > const &dxinv, Real alpha, Real beta, int ncomp) noexcept
Definition: AMReX_MLEBABecLap_2D_K.H:11
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlebabeclap_adotx(Box const &box, Array4< Real > const &y, Array4< Real const > const &x, Array4< Real const > const &a, Array4< Real const > const &bX, Array4< Real const > const &bY, Array4< const int > const &ccm, Array4< EBCellFlag const > const &flag, Array4< Real const > const &vfrc, Array4< Real const > const &apx, Array4< Real const > const &apy, Array4< Real const > const &fcx, Array4< Real const > const &fcy, Array4< Real const > const &ba, Array4< Real const > const &bc, Array4< Real const > const &beb, bool is_dirichlet, Array4< Real const > const &phieb, bool is_inhomog, GpuArray< Real, AMREX_SPACEDIM > const &dxinv, Real alpha, Real beta, int ncomp, bool beta_on_centroid, bool phi_on_centroid) noexcept
Definition: AMReX_MLEBABecLap_2D_K.H:166
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE T abs(const GpuComplex< T > &a_z) noexcept
Return the absolute value of a complex number.
Definition: AMReX_GpuComplex.H:356
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlebabeclap_flux_x_0(Box const &box, Array4< Real > const &fx, Array4< Real const > const &apx, Array4< Real const > const &sol, Array4< Real const > const &bX, Real dhx, int face_only, int ncomp, Box const &xbox) noexcept
Definition: AMReX_MLEBABecLap_2D_K.H:653
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlebabeclap_gsrb(Box const &box, Array4< Real > const &phi, Array4< Real const > const &rhs, Real alpha, Array4< Real const > const &a, Real dhx, Real dhy, Real dh, GpuArray< Real, AMREX_SPACEDIM > const &dx, Array4< Real const > const &bX, Array4< Real const > const &bY, Array4< int const > const &m0, Array4< int const > const &m2, Array4< int const > const &m1, Array4< int const > const &m3, Array4< Real const > const &f0, Array4< Real const > const &f2, Array4< Real const > const &f1, Array4< Real const > const &f3, Array4< const int > const &ccm, Array4< Real const > const &beb, EBData const &ebdata, bool is_dirichlet, bool beta_on_centroid, bool phi_on_centroid, Box const &vbox, int redblack, int ncomp) noexcept
Definition: AMReX_MLEBABecLap_2D_K.H:371
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 ubound(Array4< T > const &a) noexcept
Definition: AMReX_Array4.H:315
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlebabeclap_grad_x(Box const &box, Array4< Real > const &gx, Array4< Real const > const &sol, Array4< Real const > const &apx, Array4< Real const > const &fcx, Array4< int const > const &ccm, Real dxi, int ncomp, bool phi_on_centroid) noexcept
Definition: AMReX_MLEBABecLap_2D_K.H:691
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 lbound(Array4< T > const &a) noexcept
Definition: AMReX_Array4.H:308
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real grad_x_of_phi_on_centroids_extdir(int i, int j, int k, int n, Array4< Real const > const &phi, Array4< Real const > const &phieb, Array4< EBCellFlag const > const &flag, Array4< Real const > const &ccent, Array4< Real const > const &bcent, Array4< Real const > const &vfrac, Real &yloc_on_xface, bool is_eb_dirichlet, bool is_eb_inhomog, const bool on_x_face, const int domlo_x, const int domhi_x, const bool on_y_face, const int domlo_y, const int domhi_y)
Definition: AMReX_EB_LeastSquares_2D_K.H:495
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlebabeclap_flux_y_0(Box const &box, Array4< Real > const &fy, Array4< Real const > const &apy, Array4< Real const > const &sol, Array4< Real const > const &bY, Real dhy, int face_only, int ncomp, Box const &ybox) noexcept
Definition: AMReX_MLEBABecLap_2D_K.H:672
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real grad_eb_of_phi_on_centroids_extdir(int i, int j, int k, int n, Array4< Real const > const &phi, Array4< Real const > const &phieb, Array4< EBCellFlag const > const &flag, Array4< Real const > const &ccent, Array4< Real const > const &bcent, Array4< Real const > const &vfrac, Real &nrmx, Real &nrmy, bool is_eb_inhomog, const bool on_x_face, const int domlo_x, const int domhi_x, const bool on_y_face, const int domlo_y, const int domhi_y)
Definition: AMReX_EB_LeastSquares_2D_K.H:765
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlebabeclap_flux_x(Box const &box, Array4< Real > const &fx, Array4< Real const > const &apx, Array4< Real const > const &fcx, Array4< Real const > const &sol, Array4< Real const > const &bX, Array4< int const > const &ccm, Real dhx, int face_only, int ncomp, Box const &xbox, bool beta_on_centroid, bool phi_on_centroid) noexcept
Definition: AMReX_MLEBABecLap_2D_K.H:587
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real grad_y_of_phi_on_centroids_extdir(int i, int j, int k, int n, Array4< Real const > const &phi, Array4< Real const > const &phieb, Array4< EBCellFlag const > const &flag, Array4< Real const > const &ccent, Array4< Real const > const &bcent, Array4< Real const > const &vfrac, Real &xloc_on_yface, bool is_eb_dirichlet, bool is_eb_inhomog, const bool on_x_face, const int domlo_x, const int domhi_x, const bool on_y_face, const int domlo_y, const int domhi_y)
Definition: AMReX_EB_LeastSquares_2D_K.H:632
AMREX_GPU_HOST_DEVICE AMREX_ATTRIBUTE_FLATTEN_FOR void LoopConcurrent(Dim3 lo, Dim3 hi, F const &f) noexcept
Definition: AMReX_Loop.H:150
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlebabeclap_ebflux(int i, int j, int k, int n, Array4< Real > const &feb, Array4< Real const > const &x, Array4< EBCellFlag const > const &flag, Array4< Real const > const &vfrc, Array4< Real const > const &apx, Array4< Real const > const &apy, Array4< Real const > const &bc, Array4< Real const > const &beb, Array4< Real const > const &phieb, bool is_inhomog, GpuArray< Real, AMREX_SPACEDIM > const &dxinv) noexcept
Definition: AMReX_MLEBABecLap_2D_K.H:304
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlebabeclap_normalize(Box const &box, Array4< Real > const &phi, Real alpha, Array4< Real const > const &a, Real dhx, Real dhy, Real dh, const amrex::GpuArray< Real, AMREX_SPACEDIM > &dx, Array4< Real const > const &bX, Array4< Real const > const &bY, Array4< const int > const &ccm, Array4< EBCellFlag const > const &flag, Array4< Real const > const &vfrc, Array4< Real const > const &apx, Array4< Real const > const &apy, Array4< Real const > const &fcx, Array4< Real const > const &fcy, Array4< Real const > const &ba, Array4< Real const > const &bc, Array4< Real const > const &beb, bool is_dirichlet, bool beta_on_centroid, int ncomp) noexcept
Definition: AMReX_MLEBABecLap_2D_K.H:767
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlebabeclap_grad_y(Box const &box, Array4< Real > const &gy, Array4< Real const > const &sol, Array4< Real const > const &apy, Array4< Real const > const &fcy, Array4< int const > const &ccm, Real dyi, int ncomp, bool phi_on_centroid) noexcept
Definition: AMReX_MLEBABecLap_2D_K.H:715
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real grad_x_of_phi_on_centroids(int i, int j, int k, int n, Array4< Real const > const &phi, Array4< Real const > const &phieb, Array4< EBCellFlag const > const &flag, Array4< Real const > const &ccent, Array4< Real const > const &bcent, Real &yloc_on_xface, bool is_eb_dirichlet, bool is_eb_inhomog)
Definition: AMReX_EB_LeastSquares_2D_K.H:246
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlebabeclap_flux_y(Box const &box, Array4< Real > const &fy, Array4< Real const > const &apy, Array4< Real const > const &fcy, Array4< Real const > const &sol, Array4< Real const > const &bY, Array4< int const > const &ccm, Real dhy, int face_only, int ncomp, Box const &ybox, bool beta_on_centroid, bool phi_on_centroid) noexcept
Definition: AMReX_MLEBABecLap_2D_K.H:620
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE GpuComplex< T > sqrt(const GpuComplex< T > &a_z) noexcept
Return the square root of a complex number.
Definition: AMReX_GpuComplex.H:373
AMREX_GPU_HOST_DEVICE AMREX_ATTRIBUTE_FLATTEN_FOR void Loop(Dim3 lo, Dim3 hi, F const &f) noexcept
Definition: AMReX_Loop.H:125
Definition: AMReX_Array4.H:61
Definition: AMReX_EBData.H:26