OPAL (Object Oriented Parallel Accelerator Library) 2024.2
OPAL
BoxLibParticle.hpp
Go to the documentation of this file.
1//
2// Class BoxLibParticle
3// Particle class for AMReX. It works together with BoxLibLayout.
4// The class does the scatter and gather operations of attributes
5// to and from the grid. Ippl implements the same functionality in the
6// attribute class.
7//
8// Copyright (c) 2016 - 2020, Matthias Frey, Uldis Locans, Paul Scherrer Institut, Villigen PSI, Switzerland
9// All rights reserved
10//
11// Implemented as part of the PhD thesis
12// "Precise Simulations of Multibunches in High Intensity Cyclotrons"
13//
14// This file is part of OPAL.
15//
16// OPAL is free software: you can redistribute it and/or modify
17// it under the terms of the GNU General Public License as published by
18// the Free Software Foundation, either version 3 of the License, or
19// (at your option) any later version.
20//
21// You should have received a copy of the GNU General Public License
22// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
23//
24#ifndef BOXLIB_PARTICLE_HPP
25#define BOXLIB_PARTICLE_HPP
26
28#include "Utility/IpplTimings.h"
29
30#include <AMReX_BLFort.H>
31#include <AMReX_MultiFabUtil.H>
32#include <AMReX_MultiFabUtil_F.H>
33#include <AMReX_Interpolater.H>
34#include <AMReX_FillPatchUtil.H>
35
36template<class PLayout>
41
42
43template<class PLayout>
48
49
50template<class PLayout>
51template<class FT, unsigned Dim, class PT>
54 int lbase, int lfine,
55 const ParticleAttrib<int>& pbin, int bin)
56{
57 if ( lbase == lfine ) {
58 this->scatter(attrib, *(f[lbase].get()), pp, pbin, bin, lbase);
59 return;
60 }
61
62 const PLayout *layout_p = &this->getLayout();
63 int nGrow = layout_p->refRatio(lbase)[0];
64
65 AmrScalarFieldContainer_t tmp(lfine+1);
66 for (int lev = lbase; lev <= lfine; ++lev) {
67
68 f[lev]->setVal(0.0, f[lev]->nGrow());
69
70 const AmrGrid_t& ba = f[lev]->boxArray();
71 const AmrProcMap_t& dm = f[lev]->DistributionMap();
72 tmp[lev].reset(new AmrField_t(ba, dm, 1, nGrow));
73 tmp[lev]->setVal(0.0, nGrow);
74 }
75
76 this->AssignDensityFort(attrib, tmp, lbase, 1, lfine, pbin, bin);
77
78 for (int lev = lbase; lev <= lfine; ++lev)
79 AmrField_t::Copy(*f[lev], *tmp[lev], 0, 0, f[lev]->nComp(), f[lev]->nGrow());
80}
81
82
83template<class PLayout>
84template <class FT, unsigned Dim, class PT>
87 const ParticleAttrib<int>& pbin, int bin,
88 int level)
89{
90 const AmrGrid_t& ba = f.boxArray();
91 const AmrProcMap_t& dmap = f.DistributionMap();
92
93 AmrField_t tmp(ba, dmap, f.nComp(), 1);
94 tmp.setVal(0.0, 1);
95
96 this->AssignCellDensitySingleLevelFort(attrib, tmp, level, pbin, bin);
97
98 f.setVal(0.0, f.nGrow());
99
100 AmrField_t::Copy(f, tmp, 0, 0, f.nComp(), f.nGrow());
101}
102
103
104template<class PLayout>
105template<class FT, unsigned Dim, class PT>
108 int lbase, int lfine)
109{
110 this->InterpolateFort(attrib, f, lbase, lfine);
111}
112
113
114// Function from BoxLib adjusted to work with Ippl BoxLibParticle class
115// Scatter the particle attribute pa on the grid
116template<class PLayout>
117template <class AType>
119 AmrScalarFieldContainer_t& mf_to_be_filled,
120 int lev_min, int ncomp, int finest_level,
121 const ParticleAttrib<int>& pbin, int bin) const
122{
123// BL_PROFILE("AssignDensityFort()");
125 const PLayout *layout_p = &this->getLayout();
126
127 // not done in amrex
128 int rho_index = 0;
129
130 amrex::PhysBCFunct cphysbc, fphysbc;
131 int lo_bc[] = {INT_DIR, INT_DIR, INT_DIR}; // periodic boundaries
132 int hi_bc[] = {INT_DIR, INT_DIR, INT_DIR};
133 amrex::Vector<amrex::BCRec> bcs(1, amrex::BCRec(lo_bc, hi_bc));
134 amrex::PCInterp mapper;
135
136 AmrScalarFieldContainer_t tmp(finest_level+1);
137 for (int lev = lev_min; lev <= finest_level; ++lev) {
138 const AmrGrid_t& ba = mf_to_be_filled[lev]->boxArray();
139 const AmrProcMap_t& dm = mf_to_be_filled[lev]->DistributionMap();
140 tmp[lev].reset(new AmrField_t(ba, dm, 1, 0));
141 tmp[lev]->setVal(0.0);
142 }
143
144 for (int lev = lev_min; lev <= finest_level; ++lev) {
145 AssignCellDensitySingleLevelFort(pa, *mf_to_be_filled[lev], lev, pbin, bin, 1, 0);
146
147 if (lev < finest_level) {
148 amrex::InterpFromCoarseLevel(*tmp[lev+1], 0.0, *mf_to_be_filled[lev],
149 rho_index, rho_index, ncomp,
150 layout_p->Geom(lev), layout_p->Geom(lev+1),
151 cphysbc, fphysbc,
152 layout_p->refRatio(lev), &mapper, bcs);
153 }
154
155 if (lev > lev_min) {
156 // Note - this will double count the mass on the coarse level in
157 // regions covered by the fine level, but this will be corrected
158 // below in the call to average_down.
159 amrex::sum_fine_to_coarse(*mf_to_be_filled[lev],
160 *mf_to_be_filled[lev-1],
161 rho_index, 1, layout_p->refRatio(lev-1),
162 layout_p->Geom(lev-1), layout_p->Geom(lev));
163 }
164
165 mf_to_be_filled[lev]->plus(*tmp[lev], rho_index, ncomp, 0);
166 }
167
168 for (int lev = finest_level - 1; lev >= lev_min; --lev) {
169 amrex::average_down(*mf_to_be_filled[lev+1],
170 *mf_to_be_filled[lev], rho_index, ncomp, layout_p->refRatio(lev));
171 }
172
174}
175
176
177// This is the single-level version for cell-centered density
178template<class PLayout>
179template <class AType>
181 AmrField_t& mf_to_be_filled,
182 int level,
183 const ParticleAttrib<int>& pbin,
184 int bin,
185 int ncomp,
186 int /*particle_lvl_offset*/) const
187{
188// BL_PROFILE("ParticleContainer<NStructReal, NStructInt, NArrayReal, NArrayInt>::AssignCellDensitySingleLevelFort()");
189
190 const PLayout *layout_p = &this->getLayout();
191
192 AmrField_t* mf_pointer;
193
194 if (layout_p->OnSameGrids(level, mf_to_be_filled)) {
195 // If we are already working with the internal mf defined on the
196 // particle_box_array, then we just work with this.
197 mf_pointer = &mf_to_be_filled;
198 }
199 else {
200 // If mf_to_be_filled is not defined on the particle_box_array, then we need
201 // to make a temporary here and copy into mf_to_be_filled at the end.
202 mf_pointer = new AmrField_t(layout_p->ParticleBoxArray(level),
203 layout_p->ParticleDistributionMap(level),
204 ncomp, mf_to_be_filled.nGrow());
205 }
206
207 // We must have ghost cells for each FAB so that a particle in one grid can spread
208 // its effect to an adjacent grid by first putting the value into ghost cells of its
209 // own grid. The mf->sumBoundary call then adds the value from one grid's ghost cell
210 // to another grid's valid region.
211 if (mf_pointer->nGrow() < 1)
212 throw OpalException("BoxLibParticle::AssignCellDensitySingleLevelFort()",
213 "Must have at least one ghost cell when in AssignCellDensitySingleLevelFort");
214
215 const AmrGeometry_t& gm = layout_p->Geom(level);
216 const AmrReal_t* plo = gm.ProbLo();
217// const AmrReal_t* dx_particle = layout_p->Geom(level + particle_lvl_offset).CellSize();
218 const AmrReal_t* dx = gm.CellSize();
219
220 if (gm.isAnyPeriodic() && ! gm.isAllPeriodic()) {
221 throw OpalException("BoxLibParticle::AssignCellDensitySingleLevelFort()",
222 "Problem must be periodic in no or all directions");
223 }
224
225 for (amrex::MFIter mfi(*mf_pointer); mfi.isValid(); ++mfi) {
226 (*mf_pointer)[mfi].setVal(0);
227 }
228
229 //loop trough particles and distribute values on the grid
230 const ParticleLevelCounter_t& LocalNumPerLevel = this->getLocalNumPerLevel();
231 size_t lBegin = LocalNumPerLevel.begin(level);
232 size_t lEnd = LocalNumPerLevel.end(level);
233
234 AmrReal_t inv_dx[3] = { 1.0 / dx[0], 1.0 / dx[1], 1.0 / dx[2] };
235 double lxyz[3] = { 0.0, 0.0, 0.0 };
236 double wxyz_hi[3] = { 0.0, 0.0, 0.0 };
237 double wxyz_lo[3] = { 0.0, 0.0, 0.0 };
238 int ijk[3] = {0, 0, 0};
239
240 for (size_t ip = lBegin; ip < lEnd; ++ip) {
241
242 if ( bin > -1 && pbin[ip] != bin )
243 continue;
244
245 const int grid = this->Grid[ip];
246 FArrayBox_t& fab = (*mf_pointer)[grid];
247
248 // not callable:
249 // begin amrex_deposit_cic(pbx.data(), nstride, N, fab.dataPtr(), box.loVect(), box.hiVect(), plo, dx);
250 for (int i = 0; i < 3; ++i) {
251 lxyz[i] = ( this->R[ip](i) - plo[i] ) * inv_dx[i] + 0.5;
252 ijk[i] = lxyz[i];
253 wxyz_hi[i] = lxyz[i] - ijk[i];
254 wxyz_lo[i] = 1.0 - wxyz_hi[i];
255 }
256
257 int& i = ijk[0];
258 int& j = ijk[1];
259 int& k = ijk[2];
260
261 AmrIntVect_t i1(i-1, j-1, k-1);
262 AmrIntVect_t i2(i-1, j-1, k);
263 AmrIntVect_t i3(i-1, j, k-1);
264 AmrIntVect_t i4(i-1, j, k);
265 AmrIntVect_t i5(i, j-1, k-1);
266 AmrIntVect_t i6(i, j-1, k);
267 AmrIntVect_t i7(i, j, k-1);
268 AmrIntVect_t i8(i, j, k);
269
270 fab(i1, 0) += wxyz_lo[0]*wxyz_lo[1]*wxyz_lo[2]*pa[ip];
271 fab(i2, 0) += wxyz_lo[0]*wxyz_lo[1]*wxyz_hi[2]*pa[ip];
272 fab(i3, 0) += wxyz_lo[0]*wxyz_hi[1]*wxyz_lo[2]*pa[ip];
273 fab(i4, 0) += wxyz_lo[0]*wxyz_hi[1]*wxyz_hi[2]*pa[ip];
274 fab(i5, 0) += wxyz_hi[0]*wxyz_lo[1]*wxyz_lo[2]*pa[ip];
275 fab(i6, 0) += wxyz_hi[0]*wxyz_lo[1]*wxyz_hi[2]*pa[ip];
276 fab(i7, 0) += wxyz_hi[0]*wxyz_hi[1]*wxyz_lo[2]*pa[ip];
277 fab(i8, 0) += wxyz_hi[0]*wxyz_hi[1]*wxyz_hi[2]*pa[ip];
278 // end of amrex_deposit_cic
279 }
280
281 mf_pointer->SumBoundary(gm.periodicity());
282
283 // Only multiply the first component by (1/vol) because this converts mass
284 // to density. If there are additional components (like velocity), we don't
285 // want to divide those by volume.
286 const AmrReal_t vol = D_TERM(dx[0], *dx[1], *dx[2]);
287
288 mf_pointer->mult(1.0/vol, 0, 1, mf_pointer->nGrow());
289
290 // If mf_to_be_filled is not defined on the particle_box_array, then we need
291 // to copy here from mf_pointer into mf_to_be_filled. I believe that we don't
292 // need any information in ghost cells so we don't copy those.
293 if (mf_pointer != &mf_to_be_filled) {
294 mf_to_be_filled.copy(*mf_pointer,0,0,ncomp);
295 delete mf_pointer;
296 }
297}
298
299
300template<class PLayout>
301template <class AType>
303 AmrVectorFieldContainer_t& mesh_data,
304 int lev_min, int lev_max)
305{
306 for (int lev = lev_min; lev <= lev_max; ++lev) {
307 if ( lev > 0 )
308 InterpolateMultiLevelFort(pa, mesh_data, lev);
309 else
310 InterpolateSingleLevelFort(pa, mesh_data[lev], lev);
311 }
312}
313
314
315template<class PLayout>
316template <class AType>
318 AmrVectorField_t& mesh_data, int lev)
319{
320 for (std::size_t i = 0; i < mesh_data.size(); ++i) {
321 if (mesh_data[i]->nGrow() < 1)
322 throw OpalException("BoxLibParticle::InterpolateSingleLevelFort()",
323 "Must have at least one ghost cell when in InterpolateSingleLevelFort");
324 }
325
326 PLayout *layout_p = &this->getLayout();
327
328 const AmrGeometry_t& gm = layout_p->Geom(lev);
329 const AmrReal_t* plo = gm.ProbLo();
330 const AmrReal_t* dx = gm.CellSize();
331
332 //loop trough particles and distribute values on the grid
333 const ParticleLevelCounter_t& LocalNumPerLevel = this->getLocalNumPerLevel();
334 size_t lBegin = LocalNumPerLevel.begin(lev);
335 size_t lEnd = LocalNumPerLevel.end(lev);
336
337 // make sure that boundaries are filled!
338 for (std::size_t i = 0; i < mesh_data.size(); ++i) {
339 mesh_data[i]->FillBoundary(gm.periodicity());
340 }
341
342 AmrReal_t inv_dx[3] = { 1.0 / dx[0], 1.0 / dx[1], 1.0 / dx[2] };
343 double lxyz[3] = { 0.0, 0.0, 0.0 };
344 double wxyz_hi[3] = { 0.0, 0.0, 0.0 };
345 double wxyz_lo[3] = { 0.0, 0.0, 0.0 };
346 int ijk[3] = {0, 0, 0};
347 for (size_t ip = lBegin; ip < lEnd; ++ip) {
348
349 const int grid = this->Grid[ip];
350 FArrayBox_t& exfab = (*mesh_data[0])[grid];
351 FArrayBox_t& eyfab = (*mesh_data[1])[grid];
352 FArrayBox_t& ezfab = (*mesh_data[2])[grid];
353
354 // not callable
355 // begin amrex_interpolate_cic(pbx.data(), nstride, N, fab.dataPtr(), box.loVect(), box.hiVect(), nComp, plo, dx);
356 for (int i = 0; i < 3; ++i) {
357 lxyz[i] = ( this->R[ip](i) - plo[i] ) * inv_dx[i] + 0.5;
358 ijk[i] = lxyz[i];
359 wxyz_hi[i] = lxyz[i] - ijk[i];
360 wxyz_lo[i] = 1.0 - wxyz_hi[i];
361 }
362
363 int& i = ijk[0];
364 int& j = ijk[1];
365 int& k = ijk[2];
366
367 AmrIntVect_t i1(i-1, j-1, k-1);
368 AmrIntVect_t i2(i-1, j-1, k);
369 AmrIntVect_t i3(i-1, j, k-1);
370 AmrIntVect_t i4(i-1, j, k);
371 AmrIntVect_t i5(i, j-1, k-1);
372 AmrIntVect_t i6(i, j-1, k);
373 AmrIntVect_t i7(i, j, k-1);
374 AmrIntVect_t i8(i, j, k);
375
376 pa[ip](0) = wxyz_lo[0]*wxyz_lo[1]*wxyz_lo[2]*exfab(i1) +
377 wxyz_lo[0]*wxyz_lo[1]*wxyz_hi[2]*exfab(i2) +
378 wxyz_lo[0]*wxyz_hi[1]*wxyz_lo[2]*exfab(i3) +
379 wxyz_lo[0]*wxyz_hi[1]*wxyz_hi[2]*exfab(i4) +
380 wxyz_hi[0]*wxyz_lo[1]*wxyz_lo[2]*exfab(i5) +
381 wxyz_hi[0]*wxyz_lo[1]*wxyz_hi[2]*exfab(i6) +
382 wxyz_hi[0]*wxyz_hi[1]*wxyz_lo[2]*exfab(i7) +
383 wxyz_hi[0]*wxyz_hi[1]*wxyz_hi[2]*exfab(i8);
384
385 pa[ip](1) = wxyz_lo[0]*wxyz_lo[1]*wxyz_lo[2]*eyfab(i1) +
386 wxyz_lo[0]*wxyz_lo[1]*wxyz_hi[2]*eyfab(i2) +
387 wxyz_lo[0]*wxyz_hi[1]*wxyz_lo[2]*eyfab(i3) +
388 wxyz_lo[0]*wxyz_hi[1]*wxyz_hi[2]*eyfab(i4) +
389 wxyz_hi[0]*wxyz_lo[1]*wxyz_lo[2]*eyfab(i5) +
390 wxyz_hi[0]*wxyz_lo[1]*wxyz_hi[2]*eyfab(i6) +
391 wxyz_hi[0]*wxyz_hi[1]*wxyz_lo[2]*eyfab(i7) +
392 wxyz_hi[0]*wxyz_hi[1]*wxyz_hi[2]*eyfab(i8);
393
394 pa[ip](2) = wxyz_lo[0]*wxyz_lo[1]*wxyz_lo[2]*ezfab(i1) +
395 wxyz_lo[0]*wxyz_lo[1]*wxyz_hi[2]*ezfab(i2) +
396 wxyz_lo[0]*wxyz_hi[1]*wxyz_lo[2]*ezfab(i3) +
397 wxyz_lo[0]*wxyz_hi[1]*wxyz_hi[2]*ezfab(i4) +
398 wxyz_hi[0]*wxyz_lo[1]*wxyz_lo[2]*ezfab(i5) +
399 wxyz_hi[0]*wxyz_lo[1]*wxyz_hi[2]*ezfab(i6) +
400 wxyz_hi[0]*wxyz_hi[1]*wxyz_lo[2]*ezfab(i7) +
401 wxyz_hi[0]*wxyz_hi[1]*wxyz_hi[2]*ezfab(i8);
402 // end amrex_interpolate_cic
403 }
404}
405
406
407template<class PLayout>
408template <class AType>
410 AmrVectorFieldContainer_t& mesh_data,
411 int lev)
412{
413 for (std::size_t i = 0; i < mesh_data[lev].size(); ++i) {
414 if (mesh_data[lev][i]->nGrow() < 1)
415 throw OpalException("BoxLibParticle::InterpolateMultiLevelFort()",
416 "Must have at least one ghost cell when in InterpolateMultiLevelFort");
417 }
418
419 PLayout *layout_p = &this->getLayout();
420
421 const AmrGeometry_t& gm = layout_p->Geom(lev);
422 const AmrReal_t* plo = gm.ProbLo();
423 const AmrReal_t* fdx = gm.CellSize();
424 const AmrReal_t* cdx = layout_p->Geom(lev-1).CellSize();
425 const AmrReal_t* cplo = layout_p->Geom(lev-1).ProbLo();
426
427 //loop trough particles and distribute values on the grid
428 const ParticleLevelCounter_t& LocalNumPerLevel = this->getLocalNumPerLevel();
429 size_t lBegin = LocalNumPerLevel.begin(lev);
430 size_t lEnd = LocalNumPerLevel.end(lev);
431
432 // make sure that boundaries are filled!
433 for (std::size_t i = 0; i < mesh_data[lev].size(); ++i) {
434 mesh_data[lev][i]->FillBoundary(gm.periodicity());
435 }
436
437 AmrReal_t inv_fdx[3] = { 1.0 / fdx[0], 1.0 / fdx[1], 1.0 / fdx[2] };
438 AmrReal_t inv_cdx[3] = { 1.0 / cdx[0], 1.0 / cdx[1], 1.0 / cdx[2] };
439 double lxyz[3] = { 0.0, 0.0, 0.0 };
440 double wxyz_hi[3] = { 0.0, 0.0, 0.0 };
441 double wxyz_lo[3] = { 0.0, 0.0, 0.0 };
442 int ijk[3] = { 0, 0, 0 };
443
444 const AmrGrid_t& fba = mesh_data[lev][0]->boxArray();
445 const AmrProcMap_t& fdmap = mesh_data[lev][0]->DistributionMap();
446 AmrGrid_t cba = fba;
447 cba.coarsen(AmrIntVect_t(2, 2, 2));
448
449 AmrField_t cmesh_exdata(cba, fdmap,
450 mesh_data[lev][0]->nComp(),
451 mesh_data[lev][0]->nGrow());
452
453 cmesh_exdata.setVal(0.0, 0, 1, mesh_data[lev][0]->nGrow());
454 cmesh_exdata.copy(*mesh_data[lev-1][0], 0, 0, 1, 1, 1);
455 cmesh_exdata.FillBoundary(gm.periodicity());
456
457 AmrField_t cmesh_eydata(cba, fdmap,
458 mesh_data[lev][1]->nComp(),
459 mesh_data[lev][1]->nGrow());
460 cmesh_eydata.setVal(0.0, 0, 1, mesh_data[lev][1]->nGrow());
461 cmesh_eydata.copy(*mesh_data[lev-1][1], 0, 0, 1, 1, 1);
462 cmesh_eydata.FillBoundary(gm.periodicity());
463
464 AmrField_t cmesh_ezdata(cba, fdmap,
465 mesh_data[lev][2]->nComp(),
466 mesh_data[lev][2]->nGrow());
467 cmesh_ezdata.setVal(0.0, 0, 1, mesh_data[lev][2]->nGrow());
468 cmesh_ezdata.copy(*mesh_data[lev-1][2], 0, 0, 1, 1, 1);
469 cmesh_ezdata.FillBoundary(gm.periodicity());
470
471 for (size_t ip = lBegin; ip < lEnd; ++ip) {
472
473 const int grid = this->Grid[ip];
474
475 FArrayBox_t& exfab = (*(mesh_data[lev][0]))[grid];
476 FArrayBox_t& eyfab = (*(mesh_data[lev][1]))[grid];
477 FArrayBox_t& ezfab = (*(mesh_data[lev][2]))[grid];
478
479 FArrayBox_t& cexfab = cmesh_exdata[grid];
480 FArrayBox_t& ceyfab = cmesh_eydata[grid];
481 FArrayBox_t& cezfab = cmesh_ezdata[grid];
482
483 const typename PLayout::basefab_t& mfab = (*layout_p->getLevelMask(lev))[grid];
484
485 // not callable
486 // begin amrex_interpolate_cic(pbx.data(), nstride, N, fab.dataPtr(), box.loVect(), box.hiVect(), nComp, plo, dx);
487 for (int ii = 0; ii < 3; ++ii) {
488 lxyz[ii] = ( this->R[ip](ii) - plo[ii] ) * inv_fdx[ii] + 0.5;
489 ijk[ii] = lxyz[ii];
490 wxyz_hi[ii] = lxyz[ii] - ijk[ii];
491 wxyz_lo[ii] = 1.0 - wxyz_hi[ii];
492 }
493
494 int& i = ijk[0];
495 int& j = ijk[1];
496 int& k = ijk[2];
497
498 bool use_coarse = false;
499
500 // AMReX: electrostatic_pic_3d.f90
501 // use the coarse E near the level boundary
502 if ( mfab(AmrIntVect_t(i-1, j-1, k-1)) == 1 ) {
503 for (int ii = 0; ii < 3; ++ii) {
504 lxyz[ii] = ( this->R[ip](ii) - cplo[ii] ) * inv_cdx[ii] + 0.5;
505 ijk[ii] = lxyz[ii];
506 wxyz_hi[ii] = lxyz[ii] - ijk[ii];
507 wxyz_lo[ii] = 1.0 - wxyz_hi[ii];
508 }
509 use_coarse = true;
510 }
511
512 AmrIntVect_t i1(i-1, j-1, k-1);
513 AmrIntVect_t i3(i-1, j, k-1);
514 AmrIntVect_t i5(i, j-1, k-1);
515 AmrIntVect_t i7(i, j, k-1);
516
517 AmrIntVect_t i2(i-1, j-1, k);
518 AmrIntVect_t i4(i-1, j, k);
519 AmrIntVect_t i6(i, j-1, k);
520 AmrIntVect_t i8(i, j, k);
521 if ( use_coarse ) {
522 pa[ip](0) = wxyz_lo[0] * wxyz_lo[1] * wxyz_lo[2] * cexfab(i1) +
523 wxyz_lo[0] * wxyz_hi[1] * wxyz_lo[2] * cexfab(i3) +
524 wxyz_hi[0] * wxyz_lo[1] * wxyz_lo[2] * cexfab(i5) +
525 wxyz_hi[0] * wxyz_hi[1] * wxyz_lo[2] * cexfab(i7) +
526 wxyz_lo[0] * wxyz_lo[1] * wxyz_hi[2] * cexfab(i2) +
527 wxyz_lo[0] * wxyz_hi[1] * wxyz_hi[2] * cexfab(i4) +
528 wxyz_hi[0] * wxyz_lo[1] * wxyz_hi[2] * cexfab(i6) +
529 wxyz_hi[0] * wxyz_hi[1] * wxyz_hi[2] * cexfab(i8);
530
531 pa[ip](1) = wxyz_lo[0] * wxyz_lo[1] * wxyz_lo[2] * ceyfab(i1) +
532 wxyz_lo[0] * wxyz_hi[1] * wxyz_lo[2] * ceyfab(i3) +
533 wxyz_hi[0] * wxyz_lo[1] * wxyz_lo[2] * ceyfab(i5) +
534 wxyz_hi[0] * wxyz_hi[1] * wxyz_lo[2] * ceyfab(i7) +
535 wxyz_lo[0] * wxyz_lo[1] * wxyz_hi[2] * ceyfab(i2) +
536 wxyz_lo[0] * wxyz_hi[1] * wxyz_hi[2] * ceyfab(i4) +
537 wxyz_hi[0] * wxyz_lo[1] * wxyz_hi[2] * ceyfab(i6) +
538 wxyz_hi[0] * wxyz_hi[1] * wxyz_hi[2] * ceyfab(i8);
539
540 pa[ip](2) = wxyz_lo[0] * wxyz_lo[1] * wxyz_lo[2] * cezfab(i1) +
541 wxyz_lo[0] * wxyz_hi[1] * wxyz_lo[2] * cezfab(i3) +
542 wxyz_hi[0] * wxyz_lo[1] * wxyz_lo[2] * cezfab(i5) +
543 wxyz_hi[0] * wxyz_hi[1] * wxyz_lo[2] * cezfab(i7) +
544 wxyz_lo[0] * wxyz_lo[1] * wxyz_hi[2] * cezfab(i2) +
545 wxyz_lo[0] * wxyz_hi[1] * wxyz_hi[2] * cezfab(i4) +
546 wxyz_hi[0] * wxyz_lo[1] * wxyz_hi[2] * cezfab(i6) +
547 wxyz_hi[0] * wxyz_hi[1] * wxyz_hi[2] * cezfab(i8);
548 } else {
549 pa[ip](0) = wxyz_lo[0] * wxyz_lo[1] * wxyz_lo[2] * exfab(i1) +
550 wxyz_lo[0] * wxyz_hi[1] * wxyz_lo[2] * exfab(i3) +
551 wxyz_hi[0] * wxyz_lo[1] * wxyz_lo[2] * exfab(i5) +
552 wxyz_hi[0] * wxyz_hi[1] * wxyz_lo[2] * exfab(i7) +
553 wxyz_lo[0] * wxyz_lo[1] * wxyz_hi[2] * exfab(i2) +
554 wxyz_lo[0] * wxyz_hi[1] * wxyz_hi[2] * exfab(i4) +
555 wxyz_hi[0] * wxyz_lo[1] * wxyz_hi[2] * exfab(i6) +
556 wxyz_hi[0] * wxyz_hi[1] * wxyz_hi[2] * exfab(i8);
557
558 pa[ip](1) = wxyz_lo[0] * wxyz_lo[1] * wxyz_lo[2] * eyfab(i1) +
559 wxyz_lo[0] * wxyz_hi[1] * wxyz_lo[2] * eyfab(i3) +
560 wxyz_hi[0] * wxyz_lo[1] * wxyz_lo[2] * eyfab(i5) +
561 wxyz_hi[0] * wxyz_hi[1] * wxyz_lo[2] * eyfab(i7) +
562 wxyz_lo[0] * wxyz_lo[1] * wxyz_hi[2] * eyfab(i2) +
563 wxyz_lo[0] * wxyz_hi[1] * wxyz_hi[2] * eyfab(i4) +
564 wxyz_hi[0] * wxyz_lo[1] * wxyz_hi[2] * eyfab(i6) +
565 wxyz_hi[0] * wxyz_hi[1] * wxyz_hi[2] * eyfab(i8);
566
567 pa[ip](2) = wxyz_lo[0] * wxyz_lo[1] * wxyz_lo[2] * ezfab(i1) +
568 wxyz_lo[0] * wxyz_hi[1] * wxyz_lo[2] * ezfab(i3) +
569 wxyz_hi[0] * wxyz_lo[1] * wxyz_lo[2] * ezfab(i5) +
570 wxyz_hi[0] * wxyz_hi[1] * wxyz_lo[2] * ezfab(i7) +
571 wxyz_lo[0] * wxyz_lo[1] * wxyz_hi[2] * ezfab(i2) +
572 wxyz_lo[0] * wxyz_hi[1] * wxyz_hi[2] * ezfab(i4) +
573 wxyz_hi[0] * wxyz_lo[1] * wxyz_hi[2] * ezfab(i6) +
574 wxyz_hi[0] * wxyz_hi[1] * wxyz_hi[2] * ezfab(i8);
575 }
576 // end amrex_interpolate_cic
577 }
578}
579
580#endif
amr::AmrField_t AmrField_t
Definition PBunchDefs.h:34
bool scatter(Communicate &, InputIterator, InputIterator, RandomIterator, int *, int *, const ScatterOp &)
ParticleLayout< T, Dim > PLayout
PLayout::AmrIntVect_t AmrIntVect_t
AmrParticleBase< PLayout >::AmrField_t AmrField_t
AmrParticleBase< PLayout >::ParticleLevelCounter_t ParticleLevelCounter_t
PLayout::AmrProcMap_t AmrProcMap_t
void AssignCellDensitySingleLevelFort(ParticleAttrib< AType > &pa, AmrField_t &mf, int level, const ParticleAttrib< int > &pbin, int bin=-1, int ncomp=1, int particle_lvl_offset=0) const
PLayout::AmrGeometry_t AmrGeometry_t
IpplTimings::TimerRef AssignDensityTimer_m
void InterpolateMultiLevelFort(ParticleAttrib< AType > &pa, AmrVectorFieldContainer_t &mesh_data, int lev)
void AssignDensityFort(ParticleAttrib< AType > &pa, AmrScalarFieldContainer_t &mf_to_be_filled, int lev_min, int ncomp, int finest_level, const ParticleAttrib< int > &pbin, int bin=-1) const
AmrParticleBase< PLayout >::AmrVectorField_t AmrVectorField_t
void gather(ParticleAttrib< FT > &attrib, AmrVectorFieldContainer_t &f, ParticleAttrib< Vektor< PT, Dim > > &pp, int lbase, int lfine)
PLayout::AmrGrid_t AmrGrid_t
AmrParticleBase< PLayout >::AmrVectorFieldContainer_t AmrVectorFieldContainer_t
void InterpolateFort(ParticleAttrib< AType > &pa, AmrVectorFieldContainer_t &mesh_data, int lev_min, int lev_max)
amrex::FArrayBox FArrayBox_t
PLayout::AmrReal_t AmrReal_t
void scatter(ParticleAttrib< FT > &attrib, AmrScalarFieldContainer_t &f, ParticleAttrib< Vektor< PT, Dim > > &pp, int lbase, int lfine, const ParticleAttrib< int > &pbin, int bin=-1)
AmrParticleBase< PLayout >::AmrScalarFieldContainer_t AmrScalarFieldContainer_t
void InterpolateSingleLevelFort(ParticleAttrib< AType > &pa, AmrVectorField_t &mesh_data, int lev)
The base class for all OPAL exceptions.
ParticleIndex_t Grid
const ParticleLevelCounter_t & getLocalNumPerLevel() const
static TimerRef getTimer(const char *nm)
static void stopTimer(TimerRef t)
static void startTimer(TimerRef t)