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<EBCellFlag const> const& flag,
382  Array4<Real const> const& vfrc,
383  Array4<Real const> const& apx, Array4<Real const> const& apy,
384  Array4<Real const> const& fcx, Array4<Real const> const& fcy,
385  Array4<Real const> const& ba, Array4<Real const> const& bc,
386  Array4<Real const> const& beb,
387  bool is_dirichlet, bool beta_on_centroid, bool phi_on_centroid,
388  Box const& vbox, int redblack, int ncomp) noexcept
389 {
390  const auto vlo = amrex::lbound(vbox);
391  const auto vhi = amrex::ubound(vbox);
392 
393  amrex::Loop(box, ncomp, [=] (int i, int j, int k, int n) noexcept
394  {
395  if ((i+j+k+redblack) % 2 == 0)
396  {
397  if (flag(i,j,k).isCovered())
398  {
399  phi(i,j,k,n) = Real(0.0);
400  }
401  else
402  {
403  Real cf0 = (i == vlo.x && m0(vlo.x-1,j,k) > 0)
404  ? f0(vlo.x,j,k,n) : Real(0.0);
405  Real cf1 = (j == vlo.y && m1(i,vlo.y-1,k) > 0)
406  ? f1(i,vlo.y,k,n) : Real(0.0);
407  Real cf2 = (i == vhi.x && m2(vhi.x+1,j,k) > 0)
408  ? f2(vhi.x,j,k,n) : Real(0.0);
409  Real cf3 = (j == vhi.y && m3(i,vhi.y+1,k) > 0)
410  ? f3(i,vhi.y,k,n) : Real(0.0);
411 
412  if (flag(i,j,k).isRegular())
413  {
414  Real gamma = alpha*a(i,j,k)
415  + dhx * (bX(i+1,j,k,n) + bX(i,j,k,n))
416  + dhy * (bY(i,j+1,k,n) + bY(i,j,k,n));
417 
418  Real rho = dhx * (bX(i+1,j,k,n)*phi(i+1,j,k,n)
419  + bX(i ,j,k,n)*phi(i-1,j,k,n))
420  + dhy * (bY(i,j+1,k,n)*phi(i,j+1,k,n)
421  + bY(i,j ,k,n)*phi(i,j-1,k,n));
422 
423  Real delta = dhx*(bX(i,j,k,n)*cf0 + bX(i+1,j,k,n)*cf2)
424  + dhy*(bY(i,j,k,n)*cf1 + bY(i,j+1,k,n)*cf3);
425 
426  Real res = rhs(i,j,k,n) - (gamma*phi(i,j,k,n) - rho);
427  phi(i,j,k,n) += res/(gamma-delta);
428  }
429  else
430  {
431  Real kappa = vfrc(i,j,k);
432  Real apxm = apx(i,j,k);
433  Real apxp = apx(i+1,j,k);
434  Real apym = apy(i,j,k);
435  Real apyp = apy(i,j+1,k);
436 
437  Real fxm = -bX(i,j,k,n)*phi(i-1,j,k,n);
438  Real oxm = -bX(i,j,k,n)*cf0;
439  Real sxm = bX(i,j,k,n);
440  if (apxm != Real(0.0) && apxm != Real(1.0)) {
441  int jj = j + static_cast<int>(std::copysign(Real(1.0),fcx(i,j,k)));
442  Real fracy = (ccm(i-1,jj,k) || ccm(i,jj,k))
443  ? std::abs(fcx(i,j,k)) : Real(0.0);
444  if (!beta_on_centroid && !phi_on_centroid)
445  {
446  fxm = (Real(1.0)-fracy)*fxm +
447  fracy *bX(i,jj,k,n)*(phi(i,jj,k,n)-phi(i-1,jj,k,n));
448  }
449  else if (beta_on_centroid && !phi_on_centroid)
450  {
451  fxm = (Real(1.0)-fracy)*( -phi(i-1,j,k,n)) +
452  fracy *(phi(i,jj,k,n)-phi(i-1,jj,k,n));
453  fxm *= bX(i,j,k,n);
454  }
455  oxm = Real(0.0);
456  sxm = (Real(1.0)-fracy)*sxm;
457  }
458 
459  Real fxp = bX(i+1,j,k,n)*phi(i+1,j,k,n);
460  Real oxp = bX(i+1,j,k,n)*cf2;
461  Real sxp = -bX(i+1,j,k,n);
462  if (apxp != Real(0.0) && apxp != Real(1.0)) {
463  int jj = j + static_cast<int>(std::copysign(Real(1.0),fcx(i+1,j,k)));
464  Real fracy = (ccm(i,jj,k) || ccm(i+1,jj,k))
465  ? std::abs(fcx(i+1,j,k)) : Real(0.0);
466  if (!beta_on_centroid && !phi_on_centroid)
467  {
468  fxp = (Real(1.0)-fracy)*fxp +
469  fracy *bX(i+1,jj,k,n)*(phi(i+1,jj,k,n)-phi(i,jj,k,n));
470  }
471  else if (beta_on_centroid && !phi_on_centroid)
472  {
473  fxp = (Real(1.0)-fracy)*(phi(i+1,j,k,n) ) +
474  fracy *(phi(i+1,jj,k,n)-phi(i,jj,k,n));
475  fxp *= bX(i+1,j,k,n);
476  }
477  oxp = Real(0.0);
478  sxp = (Real(1.0)-fracy)*sxp;
479  }
480 
481  Real fym = -bY(i,j,k,n)*phi(i,j-1,k,n);
482  Real oym = -bY(i,j,k,n)*cf1;
483  Real sym = bY(i,j,k,n);
484  if (apym != Real(0.0) && apym != Real(1.0)) {
485  int ii = i + static_cast<int>(std::copysign(Real(1.0),fcy(i,j,k)));
486  Real fracx = (ccm(ii,j-1,k) || ccm(ii,j,k))
487  ? std::abs(fcy(i,j,k)) : Real(0.0);
488  if (!beta_on_centroid && !phi_on_centroid)
489  {
490  fym = (Real(1.0)-fracx)*fym +
491  fracx *bY(ii,j,k,n)*(phi(ii,j,k,n)-phi(ii,j-1,k,n));
492  }
493  else if (beta_on_centroid && !phi_on_centroid)
494  {
495  fym = (Real(1.0)-fracx)*( -phi( i,j-1,k,n)) +
496  fracx *(phi(ii,j,k,n)-phi(ii,j-1,k,n));
497  fym *= bY(i,j,k,n);
498  }
499 
500  oym = Real(0.0);
501  sym = (Real(1.0)-fracx)*sym;
502  }
503 
504  Real fyp = bY(i,j+1,k,n)*phi(i,j+1,k,n);
505  Real oyp = bY(i,j+1,k,n)*cf3;
506  Real syp = -bY(i,j+1,k,n);
507  if (apyp != Real(0.0) && apyp != Real(1.0)) {
508  int ii = i + static_cast<int>(std::copysign(Real(1.0),fcy(i,j+1,k)));
509  Real fracx = (ccm(ii,j,k) || ccm(ii,j+1,k))
510  ? std::abs(fcy(i,j+1,k)) : Real(0.0);
511  if (!beta_on_centroid && !phi_on_centroid)
512  {
513  fyp = (Real(1.0)-fracx)*fyp +
514  fracx*bY(ii,j+1,k,n)*(phi(ii,j+1,k,n)-phi(ii,j,k,n));
515  }
516  else if (beta_on_centroid && !phi_on_centroid)
517  {
518  fyp = (Real(1.0)-fracx)*(phi(i,j+1,k,n) )+
519  fracx *(phi(ii,j+1,k,n)-phi(ii,j,k,n));
520  fyp *= bY(i,j+1,k,n);
521  }
522  oyp = Real(0.0);
523  syp = (Real(1.0)-fracx)*syp;
524  }
525 
526  Real vfrcinv = (Real(1.0)/kappa);
527  Real gamma = alpha*a(i,j,k) + vfrcinv *
528  (dhx*(apxm*sxm-apxp*sxp) +
529  dhy*(apym*sym-apyp*syp));
530  Real rho = -vfrcinv *
531  (dhx*(apxm*fxm-apxp*fxp) +
532  dhy*(apym*fym-apyp*fyp));
533 
534  Real delta = -vfrcinv *
535  (dhx*(apxm*oxm-apxp*oxp) +
536  dhy*(apym*oym-apyp*oyp));
537 
538  if (is_dirichlet) {
539  Real dapx = (apxm-apxp)*dx[1];
540  Real dapy = (apym-apyp)*dx[0];
541  Real anorm = std::hypot(dapx,dapy);
542  Real anorminv = Real(1.0)/anorm;
543  Real anrmx = dapx * anorminv;
544  Real anrmy = dapy * anorminv;
545 
546  Real bctx = bc(i,j,k,0);
547  Real bcty = bc(i,j,k,1);
548  Real dx_eb = get_dx_eb(vfrc(i,j,k));
549 
550  Real dg, gx, gy, sx, sy;
551  if (std::abs(anrmx) > std::abs(anrmy)) {
552  dg = dx_eb / std::abs(anrmx);
553  } else {
554  dg = dx_eb / std::abs(anrmy);
555  }
556  gx = bctx - dg*anrmx;
557  gy = bcty - dg*anrmy;
558  sx = std::copysign(Real(1.0),anrmx);
559  sy = std::copysign(Real(1.0),anrmy);
560 
561  int ii = i - static_cast<int>(sx);
562  int jj = j - static_cast<int>(sy);
563 
564  Real phig_gamma = (Real(1.0) + gx*sx + gy*sy + gx*gy*sx*sy);
565  Real phig = ( - gx*sx - gx*gy*sx*sy) * phi(ii,j ,k,n)
566  + ( - gy*sy - gx*gy*sx*sy) * phi(i ,jj,k,n)
567  + ( + gx*gy*sx*sy) * phi(ii,jj,k,n);
568 
569  // In gsrb we are always in residual-correction form so phib = 0
570  Real dphidn = ( -phig)/dg;
571 
572  Real feb = dphidn * ba(i,j,k) * beb(i,j,k,n);
573  rho += -vfrcinv*(-dh)*feb;
574 
575  Real feb_gamma = -phig_gamma/dg * ba(i,j,k) * beb(i,j,k,n);
576  gamma += vfrcinv*(-dh)*feb_gamma;
577  }
578 
579  Real res = rhs(i,j,k,n) - (gamma*phi(i,j,k,n) - rho);
580  phi(i,j,k,n) += res/(gamma-delta);
581  }
582  }
583  }
584  });
585 }
586 
588 void mlebabeclap_flux_x (Box const& box, Array4<Real> const& fx, Array4<Real const> const& apx,
589  Array4<Real const> const& fcx, Array4<Real const> const& sol,
590  Array4<Real const> const& bX, Array4<int const> const& ccm,
591  Real dhx, int face_only, int ncomp, Box const& xbox,
592  bool beta_on_centroid, bool phi_on_centroid) noexcept
593 {
594  int lof = xbox.smallEnd(0);
595  int hif = xbox.bigEnd(0);
596  amrex::LoopConcurrent(box, ncomp, [=] (int i, int j, int k, int n) noexcept
597  {
598  if (!face_only || lof == i || hif == i) {
599  if (apx(i,j,k) == Real(0.0)) {
600  fx(i,j,k,n) = Real(0.0);
601  } else if (apx(i,j,k) == Real(1.0)) {
602  fx(i,j,k,n) = -dhx*bX(i,j,k,n)*(sol(i,j,k,n)-sol(i-1,j,k,n));
603  } else {
604  Real fxm = bX(i,j,k,n)*(sol(i,j,k,n)-sol(i-1,j,k,n));
605  int jj = j + static_cast<int>(std::copysign(Real(1.0),fcx(i,j,k)));
606  Real fracy = (ccm(i-1,jj,k) || ccm(i,jj,k)) ? std::abs(fcx(i,j,k)) : Real(0.0);
607  if (!beta_on_centroid && !phi_on_centroid)
608  {
609  fxm = (Real(1.0)-fracy)*fxm + fracy*bX(i,jj,k,n)*(sol(i,jj,k,n)-sol(i-1,jj,k,n));
610  } else if (beta_on_centroid && !phi_on_centroid) {
611  fxm = bX(i,j,k,n) * ( (Real(1.0)-fracy)*(sol(i, j,k,n)-sol(i-1, j,k,n)) +
612  fracy *(sol(i,jj,k,n)-sol(i-1,jj,k,n)) );
613  }
614  fx(i,j,k,n) = -fxm*dhx;
615  }
616  }
617  });
618 }
619 
621 void mlebabeclap_flux_y (Box const& box, Array4<Real> const& fy, Array4<Real const> const& apy,
622  Array4<Real const> const& fcy, Array4<Real const> const& sol,
623  Array4<Real const> const& bY, Array4<int const> const& ccm,
624  Real dhy, int face_only, int ncomp, Box const& ybox,
625  bool beta_on_centroid, bool phi_on_centroid) noexcept
626 {
627  int lof = ybox.smallEnd(1);
628  int hif = ybox.bigEnd(1);
629  amrex::LoopConcurrent(box, ncomp, [=] (int i, int j, int k, int n) noexcept
630  {
631  if (!face_only || lof == j || hif == j) {
632  if (apy(i,j,k) == Real(0.0)) {
633  fy(i,j,k,n) = Real(0.0);
634  } else if (apy(i,j,k) == Real(1.0)) {
635  fy(i,j,k,n) = -dhy*bY(i,j,k,n)*(sol(i,j,k,n)-sol(i,j-1,k,n));
636  } else {
637  Real fym = bY(i,j,k,n)*(sol(i,j,k,n)-sol(i,j-1,k,n));
638  int ii = i + static_cast<int>(std::copysign(Real(1.0),fcy(i,j,k)));
639  Real fracx = (ccm(ii,j-1,k) || ccm(ii,j,k)) ? std::abs(fcy(i,j,k)) : Real(0.0);
640  if (!beta_on_centroid && !phi_on_centroid)
641  {
642  fym = (Real(1.0)-fracx)*fym + fracx*bY(ii,j,k,n)*(sol(ii,j,k,n)-sol(ii,j-1,k,n));
643  } else if (beta_on_centroid && !phi_on_centroid) {
644  fym = bY(i,j,k,n) * ( (Real(1.0)-fracx)*(sol( i,j,k,n)-sol( i,j-1,k,n)) +
645  fracx *(sol(ii,j,k,n)-sol(ii,j-1,k,n)) );
646  }
647  fy(i,j,k,n) = -fym*dhy;
648  }
649  }
650  });
651 }
652 
654 void mlebabeclap_flux_x_0 (Box const& box, Array4<Real> const& fx, Array4<Real const> const& apx,
655  Array4<Real const> const& sol, Array4<Real const> const& bX,
656  Real dhx, int face_only, int ncomp, Box const& xbox) noexcept
657 {
658  int lof = xbox.smallEnd(0);
659  int hif = xbox.bigEnd(0);
660  amrex::LoopConcurrent(box, ncomp, [=] (int i, int j, int k, int n) noexcept
661  {
662  if (!face_only || lof == i || hif == i) {
663  if (apx(i,j,k) == Real(0.0)) {
664  fx(i,j,k,n) = Real(0.0);
665  } else {
666  fx(i,j,k,n) = -dhx*bX(i,j,k,n)*(sol(i,j,k,n)-sol(i-1,j,k,n));
667  }
668  }
669  });
670 }
671 
673 void mlebabeclap_flux_y_0 (Box const& box, Array4<Real> const& fy, Array4<Real const> const& apy,
674  Array4<Real const> const& sol, Array4<Real const> const& bY,
675  Real dhy, int face_only, int ncomp, Box const& ybox) noexcept
676 {
677  int lof = ybox.smallEnd(1);
678  int hif = ybox.bigEnd(1);
679  amrex::LoopConcurrent(box, ncomp, [=] (int i, int j, int k, int n) noexcept
680  {
681  if (!face_only || lof == j || hif == j) {
682  if (apy(i,j,k) == Real(0.0)) {
683  fy(i,j,k,n) = Real(0.0);
684  } else {
685  fy(i,j,k,n) = -dhy*bY(i,j,k,n)*(sol(i,j,k,n)-sol(i,j-1,k,n));
686  }
687  }
688  });
689 }
690 
692 void mlebabeclap_grad_x (Box const& box, Array4<Real> const& gx, Array4<Real const> const& sol,
693  Array4<Real const> const& apx, Array4<Real const> const& fcx,
694  Array4<int const> const& ccm,
695  Real dxi, int ncomp, bool phi_on_centroid) noexcept
696 {
697  amrex::LoopConcurrent(box, ncomp, [=] (int i, int j, int k, int n) noexcept
698  {
699  if (apx(i,j,k) == Real(0.0)) {
700  gx(i,j,k,n) = Real(0.0);
701  } else if (apx(i,j,k) == Real(1.0)) {
702  gx(i,j,k,n) = dxi*(sol(i,j,k,n)-sol(i-1,j,k,n));
703  } else {
704  Real gxm = (sol(i,j,k,n)-sol(i-1,j,k,n));
705  int jj = j + static_cast<int>(std::copysign(Real(1.0),fcx(i,j,k)));
706  Real fracy = (ccm(i-1,jj,k) || ccm(i,jj,k)) ? std::abs(fcx(i,j,k)) : Real(0.0);
707  if (!phi_on_centroid) {
708  gxm = (Real(1.0)-fracy)*gxm + fracy*(sol(i,jj,k,n)-sol(i-1,jj,k,n));
709  }
710  gx(i,j,k,n) = gxm*dxi;
711  }
712  });
713 }
714 
716 void mlebabeclap_grad_y (Box const& box, Array4<Real> const& gy, Array4<Real const> const& sol,
717  Array4<Real const> const& apy, Array4<Real const> const& fcy,
718  Array4<int const> const& ccm,
719  Real dyi, int ncomp, bool phi_on_centroid) noexcept
720 {
721  amrex::LoopConcurrent(box, ncomp, [=] (int i, int j, int k, int n) noexcept
722  {
723  if (apy(i,j,k) == Real(0.0)) {
724  gy(i,j,k,n) = Real(0.0);
725  } else if (apy(i,j,k) == Real(1.0)) {
726  gy(i,j,k,n) = dyi*(sol(i,j,k,n)-sol(i,j-1,k,n));
727  } else {
728  Real gym = (sol(i,j,k,n)-sol(i,j-1,k,n));
729  int ii = i + static_cast<int>(std::copysign(Real(1.0),fcy(i,j,k)));
730  Real fracx = (ccm(ii,j-1,k) || ccm(ii,j,k)) ? std::abs(fcy(i,j,k)) : Real(0.0);
731  if (!phi_on_centroid) {
732  gym = (Real(1.0)-fracx)*gym + fracx*(sol(ii,j,k,n)-sol(ii,j-1,k,n));
733  }
734  gy(i,j,k,n) = gym*dyi;
735  }
736  });
737 }
738 
740 void mlebabeclap_grad_x_0 (Box const& box, Array4<Real> const& gx, Array4<Real const> const& sol,
741  Array4<Real const> const& apx, Real dxi, int ncomp) noexcept
742 {
743  amrex::LoopConcurrent(box, ncomp, [=] (int i, int j, int k, int n) noexcept
744  {
745  if (apx(i,j,k) == Real(0.0)) {
746  gx(i,j,k,n) = Real(0.0);
747  } else {
748  gx(i,j,k,n) = dxi*(sol(i,j,k,n)-sol(i-1,j,k,n));
749  }
750  });
751 }
752 
754 void mlebabeclap_grad_y_0 (Box const& box, Array4<Real> const& gy, Array4<Real const> const& sol,
755  Array4<Real const> const& apy, Real dyi, int ncomp) noexcept
756 {
757  amrex::LoopConcurrent(box, ncomp, [=] (int i, int j, int k, int n) noexcept
758  {
759  if (apy(i,j,k) == Real(0.0)) {
760  gy(i,j,k,n) = Real(0.0);
761  } else {
762  gy(i,j,k,n) = dyi*(sol(i,j,k,n)-sol(i,j-1,k,n));
763  }
764  });
765 }
766 
768 void mlebabeclap_normalize (Box const& box, Array4<Real> const& phi,
769  Real alpha, Array4<Real const> const& a,
770  Real dhx, Real dhy, Real dh,
772  Array4<Real const> const& bX, Array4<Real const> const& bY,
773  Array4<const int> const& ccm, Array4<EBCellFlag const> const& flag,
774  Array4<Real const> const& vfrc,
775  Array4<Real const> const& apx, Array4<Real const> const& apy,
776  Array4<Real const> const& fcx, Array4<Real const> const& fcy,
777  Array4<Real const> const& ba, Array4<Real const> const& bc,
778  Array4<Real const> const& beb,
779  bool is_dirichlet, bool beta_on_centroid, int ncomp) noexcept
780 {
781  amrex::Loop(box, ncomp, [=] (int i, int j, int k, int n) noexcept
782  {
783  if (flag(i,j,k).isRegular())
784  {
785  phi(i,j,k,n) /= alpha*a(i,j,k) + dhx*(bX(i,j,k,n) + bX(i+1,j,k,n))
786  + dhy*(bY(i,j,k,n) + bY(i,j+1,k,n));
787  }
788  else if (flag(i,j,k).isSingleValued())
789  {
790  Real kappa = vfrc(i,j,k);
791  Real apxm = apx(i,j,k);
792  Real apxp = apx(i+1,j,k);
793  Real apym = apy(i,j,k);
794  Real apyp = apy(i,j+1,k);
795 
796  Real sxm = bX(i,j,k,n);
797  if (apxm != Real(0.0) && apxm != Real(1.0) && !beta_on_centroid) {
798  int jj = j + static_cast<int>(std::copysign(Real(1.0),fcx(i,j,k)));
799  Real fracy = (ccm(i-1,jj,k) || ccm(i,jj,k))
800  ? std::abs(fcx(i,j,k)) : Real(0.0);
801  sxm = (Real(1.0)-fracy)*sxm;
802  }
803 
804  Real sxp = -bX(i+1,j,k,n);
805  if (apxp != Real(0.0) && apxp != Real(1.0) && !beta_on_centroid) {
806  int jj = j + static_cast<int>(std::copysign(Real(1.0),fcx(i+1,j,k)));
807  Real fracy = (ccm(i,jj,k) || ccm(i+1,jj,k))
808  ? std::abs(fcx(i+1,j,k)) : Real(0.0);
809  sxp = (Real(1.0)-fracy)*sxp;
810  }
811 
812  Real sym = bY(i,j,k,n);
813  if (apym != Real(0.0) && apym != Real(1.0) && !beta_on_centroid) {
814  int ii = i + static_cast<int>(std::copysign(Real(1.0),fcy(i,j,k)));
815  Real fracx = (ccm(ii,j-1,k) || ccm(ii,j,k))
816  ? std::abs(fcy(i,j,k)) : Real(0.0);
817  sym = (Real(1.0)-fracx)*sym;
818  }
819 
820  Real syp = -bY(i,j+1,k,n);
821  if (apyp != Real(0.0) && apyp != Real(1.0) && !beta_on_centroid) {
822  int ii = i + static_cast<int>(std::copysign(Real(1.0),fcy(i,j+1,k)));
823  Real fracx = (ccm(ii,j,k) || ccm(ii,j+1,k))
824  ? std::abs(fcy(i,j+1,k)) : Real(0.0);
825  syp = (Real(1.0)-fracx)*syp;
826  }
827 
828  Real vfrcinv = (Real(1.0)/kappa);
829  Real gamma = alpha*a(i,j,k) + vfrcinv *
830  (dhx*(apxm*sxm-apxp*sxp) +
831  dhy*(apym*sym-apyp*syp));
832 
833  if (is_dirichlet) {
834  Real dapx = (apxm-apxp)*dx[1];
835  Real dapy = (apym-apyp)*dx[0];
836  Real anorm = std::hypot(dapx,dapy);
837  Real anorminv = Real(1.0)/anorm;
838  Real anrmx = dapx * anorminv;
839  Real anrmy = dapy * anorminv;
840 
841  Real bctx = bc(i,j,k,0);
842  Real bcty = bc(i,j,k,1);
843  Real dx_eb = get_dx_eb(vfrc(i,j,k));
844 
845  Real dg, gx, gy, sx, sy;
846  if (std::abs(anrmx) > std::abs(anrmy)) {
847  dg = dx_eb / std::abs(anrmx);
848  } else {
849  dg = dx_eb / std::abs(anrmy);
850  }
851  gx = bctx - dg*anrmx;
852  gy = bcty - dg*anrmy;
853  sx = std::copysign(Real(1.0),anrmx);
854  sy = std::copysign(Real(1.0),anrmy);
855 
856  Real phig_gamma = (Real(1.0) + gx*sx + gy*sy + gx*gy*sx*sy);
857  Real feb_gamma = -phig_gamma/dg * ba(i,j,k) * beb(i,j,k,n);
858  gamma += vfrcinv*(-dh)*feb_gamma;
859  }
860 
861  phi(i,j,k,n) /= gamma;
862  }
863  });
864 }
865 
866 }
867 #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:740
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:754
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:654
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:692
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:673
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:588
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< 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, 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 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:768
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:716
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:621
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