Block-Structured AMR Software Framework
 
Loading...
Searching...
No Matches
AMReX_MLNodeTensorLap_3D_K.H
Go to the documentation of this file.
1#ifndef AMREX_MLNODETENSORLAP_3D_K_H_
2#define AMREX_MLNODETENSORLAP_3D_K_H_
3#include <AMReX_Config.H>
4
5namespace amrex {
6
7namespace mlndts_detail {
8
11 int ic, int jc, int kc) noexcept
12 {
13 return (crse(ic,jc,kc)+crse(ic+1,jc,kc))*Real(0.5);
14 }
15
18 int ic, int jc, int kc) noexcept
19 {
20 return (crse(ic,jc,kc)+crse(ic,jc+1,kc))*Real(0.5);
21 }
22
25 int ic, int jc, int kc) noexcept
26 {
27 return (crse(ic,jc,kc)+crse(ic,jc,kc+1))*Real(0.5);
28 }
29
32 int ic, int jc, int kc) noexcept
33 {
34 return (ts_interp_line_y(crse,ic ,jc ,kc) +
35 ts_interp_line_y(crse,ic+1,jc ,kc) +
36 ts_interp_line_x(crse,ic ,jc ,kc) +
37 ts_interp_line_x(crse,ic ,jc+1,kc)) * Real(0.25);
38 }
39
42 int ic, int jc, int kc) noexcept
43 {
44 return (ts_interp_line_z(crse,ic ,jc,kc ) +
45 ts_interp_line_z(crse,ic+1,jc,kc ) +
46 ts_interp_line_x(crse,ic ,jc,kc ) +
47 ts_interp_line_x(crse,ic ,jc,kc+1)) * Real(0.25);
48 }
49
52 int ic, int jc, int kc) noexcept
53 {
54 return (ts_interp_line_z(crse,ic,jc ,kc ) +
55 ts_interp_line_z(crse,ic,jc+1,kc ) +
56 ts_interp_line_y(crse,ic,jc ,kc ) +
57 ts_interp_line_y(crse,ic,jc ,kc+1)) * Real(0.25);
58 }
59}
60
62void mlndtslap_interpadd (int i, int j, int k, Array4<Real> const& fine,
63 Array4<Real const> const& crse, Array4<int const> const& msk) noexcept
64{
65 using namespace mlndts_detail;
66
67 if (!msk(i,j,k)) {
68 int ic = amrex::coarsen(i,2);
69 int jc = amrex::coarsen(j,2);
70 int kc = amrex::coarsen(k,2);
71 bool i_is_odd = (ic*2 != i);
72 bool j_is_odd = (jc*2 != j);
73 bool k_is_odd = (kc*2 != k);
74 if (i_is_odd && j_is_odd && k_is_odd) {
75 // Fine node at center of cell
76 fine(i,j,k) += (ts_interp_face_yz(crse,ic ,jc ,kc ) +
77 ts_interp_face_yz(crse,ic+1,jc ,kc ) +
78 ts_interp_face_xz(crse,ic ,jc ,kc ) +
79 ts_interp_face_xz(crse,ic ,jc+1,kc ) +
80 ts_interp_face_xy(crse,ic ,jc ,kc ) +
81 ts_interp_face_xy(crse,ic ,jc ,kc+1)) * Real(1./6.);
82 } else if (j_is_odd && k_is_odd) {
83 // Node on a Y-Z face
84 fine(i,j,k) += ts_interp_face_yz(crse,ic,jc,kc);
85 } else if (i_is_odd && k_is_odd) {
86 // Node on a Z-X face
87 fine(i,j,k) += ts_interp_face_xz(crse,ic,jc,kc);
88 } else if (i_is_odd && j_is_odd) {
89 // Node on a X-Y face
90 fine(i,j,k) += ts_interp_face_xy(crse,ic,jc,kc);
91 } else if (i_is_odd) {
92 // Node on X line
93 fine(i,j,k) += ts_interp_line_x(crse,ic,jc,kc);
94 } else if (j_is_odd) {
95 // Node on Y line
96 fine(i,j,k) += ts_interp_line_y(crse,ic,jc,kc);
97 } else if (k_is_odd) {
98 // Node on Z line
99 fine(i,j,k) += ts_interp_line_z(crse,ic,jc,kc);
100 } else {
101 // Node coincident with coarse node
102 fine(i,j,k) += crse(ic,jc,kc);
103 }
104 }
105}
106
108void mlndtslap_semi_interpadd (int i, int j, int k, Array4<Real> const& fine,
109 Array4<Real const> const& crse, Array4<int const> const& msk,
110 int semi_dir) noexcept
111{
112 using namespace mlndts_detail;
113
114 if (!msk(i,j,k)) {
115 if (semi_dir == 0) {
116 int jc = amrex::coarsen(j,2);
117 int kc = amrex::coarsen(k,2);
118 bool j_is_odd = (jc*2 != j);
119 bool k_is_odd = (kc*2 != k);
120 if (j_is_odd && k_is_odd) {
121 // Node on a Y-Z face
122 fine(i,j,k) += ts_interp_face_yz(crse,i,jc,kc);
123 } else if (j_is_odd) {
124 // Node on Y line
125 fine(i,j,k) += ts_interp_line_y(crse,i,jc,kc);
126 } else if (k_is_odd) {
127 // Node on Z line
128 fine(i,j,k) += ts_interp_line_z(crse,i,jc,kc);
129 } else {
130 // Node coincident with coarse node
131 fine(i,j,k) += crse(i,jc,kc);
132 }
133 } else if (semi_dir == 1) {
134 int ic = amrex::coarsen(i,2);
135 int kc = amrex::coarsen(k,2);
136 bool i_is_odd = (ic*2 != i);
137 bool k_is_odd = (kc*2 != k);
138 if (i_is_odd && k_is_odd) {
139 // Node on a Z-X face
140 fine(i,j,k) += ts_interp_face_xz(crse,ic,j,kc);
141 } else if (i_is_odd) {
142 // Node on X line
143 fine(i,j,k) += ts_interp_line_x(crse,ic,j,kc);
144 } else if (k_is_odd) {
145 // Node on Z line
146 fine(i,j,k) += ts_interp_line_z(crse,ic,j,kc);
147 } else {
148 // Node coincident with coarse node
149 fine(i,j,k) += crse(ic,j,kc);
150 }
151 } else {
152 int ic = amrex::coarsen(i,2);
153 int jc = amrex::coarsen(j,2);
154 bool i_is_odd = (ic*2 != i);
155 bool j_is_odd = (jc*2 != j);
156 if (i_is_odd && j_is_odd) {
157 // Node on a X-Y face
158 fine(i,j,k) += ts_interp_face_xy(crse,ic,jc,k);
159 } else if (i_is_odd) {
160 // Node on X line
161 fine(i,j,k) += ts_interp_line_x(crse,ic,jc,k);
162 } else if (j_is_odd) {
163 // Node on Y line
164 fine(i,j,k) += ts_interp_line_y(crse,ic,jc,k);
165 } else {
166 // Node coincident with coarse node
167 fine(i,j,k) += crse(ic,jc,k);
168 }
169 }
170 }
171}
172
174void mlndtslap_adotx (int i, int j, int k, Array4<Real> const& y, Array4<Real const> const& x,
175 Array4<int const> const& msk, GpuArray<Real,6> const& s) noexcept
176{
177 if (msk(i,j,k)) {
178 y(i,j,k) = Real(0.0);
179 } else {
180 y(i,j,k) = s[0] * (x(i-1,j ,k ) + x(i+1,j ,k ))
181 + s[3] * (x(i ,j-1,k ) + x(i ,j+1,k ))
182 + s[5] * (x(i ,j ,k-1) + x(i ,j ,k+1))
183 - Real(2.)*(s[0]+s[3]+s[5]) * x(i,j,k)
184 + Real(0.5)*s[1] * (x(i-1,j-1,k ) + x(i+1,j+1,k ) - x(i-1,j+1,k ) - x(i+1,j-1,k ))
185 + Real(0.5)*s[2] * (x(i-1,j ,k-1) + x(i+1,j ,k+1) - x(i-1,j ,k+1) - x(i+1,j ,k-1))
186 + Real(0.5)*s[4] * (x(i ,j-1,k-1) + x(i ,j+1,k+1) - x(i ,j-1,k+1) - x(i ,j+1,k-1));
187 }
188}
189
191void mlndtslap_gauss_seidel (int i, int j, int k, Array4<Real> const& sol,
192 Array4<Real const> const& rhs, Array4<int const> const& msk,
193 GpuArray<Real,6> const& s) noexcept
194{
195 if (msk(i,j,k)) {
196 sol(i,j,k) = 0.0;
197 } else {
198 constexpr Real omega = Real(1.25);
199 Real s0 = Real(-2.)*(s[0]+s[3]+s[5]);
200 Real Ax = s[0] * (sol(i-1,j ,k ) + sol(i+1,j ,k ))
201 + s[3] * (sol(i ,j-1,k ) + sol(i ,j+1,k ))
202 + s[5] * (sol(i ,j ,k-1) + sol(i ,j ,k+1))
203 + s0 * sol(i,j,k)
204 + Real(0.5)*s[1] * (sol(i-1,j-1,k ) + sol(i+1,j+1,k ) - sol(i-1,j+1,k ) - sol(i+1,j-1,k ))
205 + Real(0.5)*s[2] * (sol(i-1,j ,k-1) + sol(i+1,j ,k+1) - sol(i-1,j ,k+1) - sol(i+1,j ,k-1))
206 + Real(0.5)*s[4] * (sol(i ,j-1,k-1) + sol(i ,j+1,k+1) - sol(i ,j-1,k+1) - sol(i ,j+1,k-1));
207 sol(i,j,k) += (rhs(i,j,k) - Ax) * (omega/s0);
208 }
209}
210
211#if defined(AMREX_USE_HYPRE) && (AMREX_SPACEDIM > 1)
212
213template <typename HypreInt, typename AtomicInt>
214void mlndtslap_fill_ijmatrix_cpu (Box const& ndbx,
215 Array4<AtomicInt const> const& gid,
216 Array4<int const> const& lid,
217 HypreInt* const ncols, HypreInt* const cols,
218 Real* const mat, // NOLINT(readability-non-const-parameter)
219 GpuArray<Real,6> const& s) noexcept
220{
221 constexpr auto gidmax = std::numeric_limits<AtomicInt>::max();
222 HypreInt nelems = 0;
223 amrex::LoopOnCpu(ndbx, [&] (int i, int j, int k) noexcept
224 {
225 if (lid(i,j,k) >= 0)
226 {
227 HypreInt nelems_old = nelems;
228
229 cols[nelems] = gid(i,j,k);
230 mat[nelems] = Real(-2.)*(s[0]+s[3]+s[5]);
231 ++nelems;
232
233 if ( gid(i ,j-1,k-1) < gidmax) {
234 cols[nelems] = gid(i ,j-1,k-1);
235 mat[nelems] = Real(0.5)*s[4];
236 ++nelems;
237 }
238
239 if ( gid(i-1,j ,k-1) < gidmax) {
240 cols[nelems] = gid(i-1,j ,k-1);
241 mat[nelems] = Real(0.5)*s[2];
242 ++nelems;
243 }
244
245 if ( gid(i ,j ,k-1) < gidmax) {
246 cols[nelems] = gid(i ,j ,k-1);
247 mat[nelems] = s[5];
248 ++nelems;
249 }
250
251 if ( gid(i+1,j ,k-1) < gidmax) {
252 cols[nelems] = gid(i+1,j ,k-1);
253 mat[nelems] = Real(-0.5)*s[2];
254 ++nelems;
255 }
256
257 if ( gid(i ,j+1,k-1) < gidmax) {
258 cols[nelems] = gid(i ,j+1,k-1);
259 mat[nelems] = Real(-0.5)*s[4];
260 ++nelems;
261 }
262
263 if ( gid(i-1,j-1,k ) < gidmax) {
264 cols[nelems] = gid(i-1,j-1,k );
265 mat[nelems] = Real(0.5)*s[1];
266 ++nelems;
267 }
268
269 if ( gid(i ,j-1,k ) < gidmax) {
270 cols[nelems] = gid(i ,j-1,k );
271 mat[nelems] = s[3];
272 ++nelems;
273 }
274
275 if ( gid(i+1,j-1,k ) < gidmax) {
276 cols[nelems] = gid(i+1,j-1,k );
277 mat[nelems] = Real(-0.5)*s[1];
278 ++nelems;
279 }
280
281 if ( gid(i-1,j ,k ) < gidmax) {
282 cols[nelems] = gid(i-1,j ,k );
283 mat[nelems] = s[0];
284 ++nelems;
285 }
286
287 if ( gid(i+1,j ,k ) < gidmax) {
288 cols[nelems] = gid(i+1,j ,k );
289 mat[nelems] = s[0];
290 ++nelems;
291 }
292
293 if ( gid(i-1,j+1,k ) < gidmax) {
294 cols[nelems] = gid(i-1,j+1,k );
295 mat[nelems] = Real(-0.5)*s[1];
296 ++nelems;
297 }
298
299 if ( gid(i ,j+1,k ) < gidmax) {
300 cols[nelems] = gid(i ,j+1,k );
301 mat[nelems] = s[3];
302 ++nelems;
303 }
304
305 if ( gid(i+1,j+1,k ) < gidmax) {
306 cols[nelems] = gid(i+1,j+1,k );
307 mat[nelems] = Real(0.5)*s[1];
308 ++nelems;
309 }
310
311 if ( gid(i ,j-1,k+1) < gidmax) {
312 cols[nelems] = gid(i ,j-1,k+1);
313 mat[nelems] = Real(-0.5)*s[4];
314 ++nelems;
315 }
316
317 if ( gid(i-1,j ,k+1) < gidmax) {
318 cols[nelems] = gid(i-1,j ,k+1);
319 mat[nelems] = Real(-0.5)*s[2];
320 ++nelems;
321 }
322
323 if ( gid(i ,j ,k+1) < gidmax) {
324 cols[nelems] = gid(i ,j ,k+1);
325 mat[nelems] = s[5];
326 ++nelems;
327 }
328
329 if ( gid(i+1,j ,k+1) < gidmax) {
330 cols[nelems] = gid(i+1,j ,k+1);
331 mat[nelems] = Real(0.5)*s[2];
332 ++nelems;
333 }
334
335 if ( gid(i ,j+1,k+1) < gidmax) {
336 cols[nelems] = gid(i ,j+1,k+1);
337 mat[nelems] = Real(0.5)*s[4];
338 ++nelems;
339 }
340
341 ncols[lid(i,j,k)] = nelems - nelems_old;
342 }
343 });
344}
345
346#ifdef AMREX_USE_GPU
347
348template <typename HypreInt, typename AtomicInt>
350void mlndtslap_fill_ijmatrix_gpu (const int ps, const int i, const int j, const int k,
351 const int offset, Box const& ndbx,
352 Array4<AtomicInt const> const& gid,
353 Array4<int const> const& lid,
354 HypreInt* const ncols, HypreInt* const cols,
355 Real* const mat, // NOLINT(readability-non-const-parameter)
356 GpuArray<Real,6> const& s) noexcept
357{
358 if (lid(i,j,k) >= 0)
359 {
360 constexpr auto gidmax = std::numeric_limits<AtomicInt>::max();
361
362 if (offset == 0) {
363 cols[ps] = gid(i,j,k);
364 mat[ps] = Real(-2.)*(s[0]+s[3]+s[5]);
365 int nc = 1;
366 if (gid(i-1,j-1,k-1) < gidmax) { ++nc; }
367 if (gid(i ,j-1,k-1) < gidmax) { ++nc; }
368 if (gid(i+1,j-1,k-1) < gidmax) { ++nc; }
369 if (gid(i-1,j ,k-1) < gidmax) { ++nc; }
370 if (gid(i ,j ,k-1) < gidmax) { ++nc; }
371 if (gid(i+1,j ,k-1) < gidmax) { ++nc; }
372 if (gid(i-1,j+1,k-1) < gidmax) { ++nc; }
373 if (gid(i ,j+1,k-1) < gidmax) { ++nc; }
374 if (gid(i+1,j+1,k-1) < gidmax) { ++nc; }
375 if (gid(i-1,j-1,k ) < gidmax) { ++nc; }
376 if (gid(i ,j-1,k ) < gidmax) { ++nc; }
377 if (gid(i+1,j-1,k ) < gidmax) { ++nc; }
378 if (gid(i-1,j ,k ) < gidmax) { ++nc; }
379 if (gid(i+1,j ,k ) < gidmax) { ++nc; }
380 if (gid(i-1,j+1,k ) < gidmax) { ++nc; }
381 if (gid(i ,j+1,k ) < gidmax) { ++nc; }
382 if (gid(i+1,j+1,k ) < gidmax) { ++nc; }
383 if (gid(i-1,j-1,k+1) < gidmax) { ++nc; }
384 if (gid(i ,j-1,k+1) < gidmax) { ++nc; }
385 if (gid(i+1,j-1,k+1) < gidmax) { ++nc; }
386 if (gid(i-1,j ,k+1) < gidmax) { ++nc; }
387 if (gid(i ,j ,k+1) < gidmax) { ++nc; }
388 if (gid(i+1,j ,k+1) < gidmax) { ++nc; }
389 if (gid(i-1,j+1,k+1) < gidmax) { ++nc; }
390 if (gid(i ,j+1,k+1) < gidmax) { ++nc; }
391 if (gid(i+1,j+1,k+1) < gidmax) { ++nc; }
392 ncols[lid(i,j,k)] = nc;
393 }
394 else if (offset == 1 && gid(i-1,j-1,k-1) < gidmax) {
395 cols[ps] = gid(i-1,j-1,k-1);
396 mat[ps] = Real(0.0);
397 }
398 else if (offset == 2 && gid(i ,j-1,k-1) < gidmax) {
399 cols[ps] = gid(i ,j-1,k-1);
400 mat[ps] = Real(0.5)*s[4];
401 }
402 else if (offset == 3 && gid(i+1,j-1,k-1) < gidmax) {
403 cols[ps] = gid(i+1,j-1,k-1);
404 mat[ps] = Real(0.0);
405 }
406 else if (offset == 4 && gid(i-1,j ,k-1) < gidmax) {
407 cols[ps] = gid(i-1,j ,k-1);
408 mat[ps] = Real(0.5)*s[2];
409 }
410 else if (offset == 5 && gid(i ,j ,k-1) < gidmax) {
411 cols[ps] = gid(i ,j ,k-1);
412 mat[ps] = s[5];
413 }
414 else if (offset == 6 && gid(i+1,j ,k-1) < gidmax) {
415 cols[ps] = gid(i+1,j ,k-1);
416 mat[ps] = Real(-0.5)*s[2];
417 }
418 else if (offset == 7 && gid(i-1,j+1,k-1) < gidmax) {
419 cols[ps] = gid(i-1,j+1,k-1);
420 mat[ps] = Real(0.0);
421 }
422 else if (offset == 8 && gid(i ,j+1,k-1) < gidmax) {
423 cols[ps] = gid(i ,j+1,k-1);
424 mat[ps] = Real(-0.5)*s[4];
425 }
426 else if (offset == 9 && gid(i+1,j+1,k-1) < gidmax) {
427 cols[ps] = gid(i+1,j+1,k-1);
428 mat[ps] = Real(0.0);
429 }
430 else if (offset == 10 && gid(i-1,j-1,k ) < gidmax) {
431 cols[ps] = gid(i-1,j-1,k );
432 mat[ps] = Real(0.5)*s[1];
433 }
434 else if (offset == 11 && gid(i ,j-1,k ) < gidmax) {
435 cols[ps] = gid(i ,j-1,k );
436 mat[ps] = s[3];
437 }
438 else if (offset == 12 && gid(i+1,j-1,k ) < gidmax) {
439 cols[ps] = gid(i+1,j-1,k );
440 mat[ps] = Real(-0.5)*s[1];
441 }
442 else if (offset == 13 && gid(i-1,j ,k ) < gidmax) {
443 cols[ps] = gid(i-1,j ,k );
444 mat[ps] = s[0];
445 }
446 else if (offset == 14 && gid(i+1,j ,k ) < gidmax) {
447 cols[ps] = gid(i+1,j ,k );
448 mat[ps] = s[0];
449 }
450 else if (offset == 15 && gid(i-1,j+1,k ) < gidmax) {
451 cols[ps] = gid(i-1,j+1,k );
452 mat[ps] = Real(-0.5)*s[1];
453 }
454 else if (offset == 16 && gid(i ,j+1,k ) < gidmax) {
455 cols[ps] = gid(i ,j+1,k );
456 mat[ps] = s[3];
457 }
458 else if (offset == 17 && gid(i+1,j+1,k ) < gidmax) {
459 cols[ps] = gid(i+1,j+1,k );
460 mat[ps] = Real(0.5)*s[1];
461 }
462 else if (offset == 18 && gid(i-1,j-1,k+1) < gidmax) {
463 cols[ps] = gid(i-1,j-1,k+1);
464 mat[ps] = Real(0.0);
465 }
466 else if (offset == 19 && gid(i ,j-1,k+1) < gidmax) {
467 cols[ps] = gid(i ,j-1,k+1);
468 mat[ps] = Real(-0.5)*s[4];
469 }
470 else if (offset == 20 && gid(i+1,j-1,k+1) < gidmax) {
471 cols[ps] = gid(i+1,j-1,k+1);
472 mat[ps] = Real(0.0);
473 }
474 else if (offset == 21 && gid(i-1,j ,k+1) < gidmax) {
475 cols[ps] = gid(i-1,j ,k+1);
476 mat[ps] = Real(-0.5)*s[2];
477 }
478 else if (offset == 22 && gid(i ,j ,k+1) < gidmax) {
479 cols[ps] = gid(i ,j ,k+1);
480 mat[ps] = s[5];
481 }
482 else if (offset == 23 && gid(i+1,j ,k+1) < gidmax) {
483 cols[ps] = gid(i+1,j ,k+1);
484 mat[ps] = Real(0.5)*s[2];
485 }
486 else if (offset == 24 && gid(i-1,j+1,k+1) < gidmax) {
487 cols[ps] = gid(i-1,j+1,k+1);
488 mat[ps] = Real(0.0);
489 }
490 else if (offset == 25 && gid(i ,j+1,k+1) < gidmax) {
491 cols[ps] = gid(i ,j+1,k+1);
492 mat[ps] = Real(0.5)*s[4];
493 }
494 else if (offset == 26 && gid(i+1,j+1,k+1) < gidmax) {
495 cols[ps] = gid(i+1,j+1,k+1);
496 mat[ps] = Real(0.0);
497 }
498 }
499}
500
501#endif
502
503#endif
504
505}
506
507#endif
#define AMREX_FORCE_INLINE
Definition AMReX_Extension.H:119
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
#define AMREX_GPU_HOST_DEVICE
Definition AMReX_GpuQualifiers.H:20
Array4< int const > offset
Definition AMReX_HypreMLABecLap.cpp:1089
Array4< Real > fine
Definition AMReX_InterpFaceRegister.cpp:90
Array4< Real const > crse
Definition AMReX_InterpFaceRegister.cpp:92
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real ts_interp_line_x(Array4< Real const > const &crse, int ic, int jc) noexcept
Definition AMReX_MLNodeTensorLap_2D_K.H:10
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real ts_interp_line_z(Array4< Real const > const &crse, int ic, int jc, int kc) noexcept
Definition AMReX_MLNodeTensorLap_3D_K.H:24
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real ts_interp_line_y(Array4< Real const > const &crse, int ic, int jc) noexcept
Definition AMReX_MLNodeTensorLap_2D_K.H:16
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real ts_interp_face_xz(Array4< Real const > const &crse, int ic, int jc, int kc) noexcept
Definition AMReX_MLNodeTensorLap_3D_K.H:41
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real ts_interp_face_xy(Array4< Real const > const &crse, int ic, int jc) noexcept
Definition AMReX_MLNodeTensorLap_2D_K.H:22
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real ts_interp_face_yz(Array4< Real const > const &crse, int ic, int jc, int kc) noexcept
Definition AMReX_MLNodeTensorLap_3D_K.H:51
Definition AMReX_Amr.cpp:49
AMREX_ATTRIBUTE_FLATTEN_FOR void LoopOnCpu(Dim3 lo, Dim3 hi, F const &f) noexcept
Definition AMReX_Loop.H:355
BoxND< AMREX_SPACEDIM > Box
Definition AMReX_BaseFwd.H:27
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlndtslap_adotx(int i, int j, int k, Array4< Real > const &y, Array4< Real const > const &x, Array4< int const > const &msk, GpuArray< Real, 3 > const &s) noexcept
Definition AMReX_MLNodeTensorLap_2D_K.H:91
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlndtslap_gauss_seidel(int i, int j, int k, Array4< Real > const &sol, Array4< Real const > const &rhs, Array4< int const > const &msk, GpuArray< Real, 3 > const &s) noexcept
Definition AMReX_MLNodeTensorLap_2D_K.H:105
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > coarsen(const BoxND< dim > &b, int ref_ratio) noexcept
Coarsen BoxND by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo/rati...
Definition AMReX_Box.H:1304
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlndtslap_semi_interpadd(int, int, int, Array4< Real > const &, Array4< Real const > const &, Array4< int const > const &, int) noexcept
Definition AMReX_MLNodeTensorLap_1D_K.H:13
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mlndtslap_interpadd(int, int, int, Array4< Real > const &, Array4< Real const > const &, Array4< int const > const &) noexcept
Definition AMReX_MLNodeTensorLap_1D_K.H:8
Definition AMReX_Array4.H:61
Definition AMReX_Array.H:34