1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
|
# Copyright 2017-2023 David Grote
#
# This file is part of WarpX.
#
# License: BSD-3-Clause-LBNL
"""Provides wrappers around multiFABs
Available routines:
ExWrapper, EyWrapper, EzWrapper
BxWrapper, ByWrapper, BzWrapper
JxWrapper, JyWrapper, JzWrapper
ExFPWrapper, EyFPWrapper, EzFPWrapper
BxFPWrapper, ByFPWrapper, BzFPWrapper
JxFPWrapper, JyFPWrapper, JzFPWrapper
RhoFPWrapper, PhiFPWrapper
FFPWrapper, GFPWrapper
AxFPWrapper, AyFPWrapper, AzFPWrapper
ExCPWrapper, EyCPWrapper, EzCPWrapper
BxCPWrapper, ByCPWrapper, BzCPWrapper
JxCPWrapper, JyCPWrapper, JzCPWrapper
RhoCPWrapper
FCPWrapper, GCPWrapper
EdgeLengthsxWrapper, EdgeLengthsyWrapper, EdgeLengthszWrapper
FaceAreasxWrapper, FaceAreasyWrapper, FaceAreaszWrapper
ExFPPMLWrapper, EyFPPMLWrapper, EzFPPMLWrapper
BxFPPMLWrapper, ByFPPMLWrapper, BzFPPMLWrapper
JxFPPMLWrapper, JyFPPMLWrapper, JzFPPMLWrapper
JxFPAmpereWrapper, JyFPAmpereWrapper, JzFPAmpereWrapper
FFPPMLWrapper, GFPPMLWrapper
ExCPPMLWrapper, EyCPPMLWrapper, EzCPPMLWrapper
BxCPPMLWrapper, ByCPPMLWrapper, BzCPPMLWrapper
JxCPPMLWrapper, JyCPPMLWrapper, JzCPPMLWrapper
FCPPMLWrapper, GCPPMLWrapper
"""
import numpy as np
try:
import cupy as cp
except ImportError:
cp = None
try:
from mpi4py import MPI as mpi
comm_world = mpi.COMM_WORLD
npes = comm_world.Get_size()
except ImportError:
npes = 1
from ._libwarpx import libwarpx
class _MultiFABWrapper(object):
"""Wrapper around MultiFabs
This provides a convenient way to query and set data in the MultiFabs.
The indexing is based on global indices.
Parameters
----------
mf_name: string
The name of the MultiFab to be accessed, the tag specified when the
MultiFab is allocated
level: int
The refinement level
include_ghosts: bool, default=False
Whether to include the ghost cells.
Note that when True, the first n-ghost negative indices will refer to the lower
ghost cells.
"""
def __init__(self, mf_name, level, include_ghosts=False):
self.mf_name = mf_name
self.level = level
self.include_ghosts = include_ghosts
self.dim = libwarpx.dim
# The overlaps list is one along the axes where the grid boundaries overlap the neighboring grid,
# which is the case with node centering.
ix_type = self.mf.box_array().ix_type()
self.overlaps = self._get_indices([int(ix_type.node_centered(i)) for i in range(self.dim)], 0)
def __len__(self):
"Returns the number of blocks"
return self.mf.size
def __iter__(self):
"The iteration is over the MultiFab"
return self.mf.__iter__()
@property
def mf(self):
# Always fetch this anew in case the C++ MultiFab is recreated
warpx = libwarpx.libwarpx_so.get_instance()
# All MultiFab names have the level suffix
return warpx.multifab(f'{self.mf_name}[level={self.level}]')
@property
def shape(self):
"""Returns the shape of the global array
"""
min_box = self.mf.box_array().minimal_box()
shape = list(min_box.size - min_box.small_end)
if self.include_ghosts:
nghosts = self.mf.n_grow_vect()
shape = [shape[i] + 2*nghosts[i] for i in range(self.dim)]
shape.append(self.mf.nComp)
return tuple(shape)
def mesh(self, direction):
"""Returns the mesh along the specified direction with the appropriate centering.
Parameters
----------
direction: string
In 3d, one of 'x', 'y', or 'z'.
In 2d, Cartesian, one of 'x', or 'z'.
In RZ, one of 'r', or 'z'
In Z, 'z'.
"""
try:
if libwarpx.geometry_dim == '3d':
idir = ['x', 'y', 'z'].index(direction)
celldir = idir
elif libwarpx.geometry_dim == '2d':
idir = ['x', 'z'].index(direction)
celldir = 2*idir
elif libwarpx.geometry_dim == 'rz':
idir = ['r', 'z'].index(direction)
celldir = 2*idir
elif libwarpx.geometry_dim == '1d':
idir = ['z'].index(direction)
celldir = idir
except ValueError:
raise Exception('Inappropriate direction given')
min_box = self.mf.box_array().minimal_box()
ilo = min_box.small_end[idir]
ihi = min_box.big_end[idir]
if self.include_ghosts:
# The ghost cells are added to the upper and lower end of the global domain.
nghosts = self.mf.n_grow_vect()
ilo = list(ilo)
ihi = list(ihi)
min_box = self.mf.box_array().minimal_box()
imax = min_box.big_end
for i in range(self.dim):
if ilo[i] == 0:
ilo[i] -= nghosts[i]
if ihi[i] == imax[i]:
ihi[i] += nghosts[i]
# Cell size in the direction
warpx = libwarpx.libwarpx_so.get_instance()
dd = warpx.Geom(self.level).data().CellSize(idir)
# The centering shift
ix_type = self.mf.box_array().ix_type()
if ix_type.node_centered(idir):
# node centered
shift = 0.
else:
# cell centered
shift = 0.5*dd
lo = warpx.Geom(self.level).ProbLo(idir)
return lo + np.arange(ilo,ihi+1)*dd + shift
def _get_indices(self, index, missing):
"""Expand the index list to length three.
Parameters
----------
index: sequence of length dims
The indices for each dim
missing:
The value used to fill in the extra dimensions added
"""
result = []
for i in range(self.dim):
result.append(index[i])
for i in range(self.dim, 3):
result.append(missing)
return result
def _get_n_ghosts(self):
"""Return the list of number of ghosts. This includes the component dimension."""
nghosts = list(self._get_indices(self.mf.n_grow_vect(), 0))
# The components always has nghosts = 0
nghosts.append(0)
return nghosts
def _get_min_indices(self):
"""Returns the minimum indices, expanded to length 3"""
min_box = self.mf.box_array().minimal_box()
if self.include_ghosts:
min_box.grow(self.mf.n_grow_vect())
imin = self._get_indices(min_box.small_end, 0)
return imin
def _get_max_indices(self):
"""Returns the maximum indices, expanded to length 3.
"""
min_box = self.mf.box_array().minimal_box()
if self.include_ghosts:
min_box.grow(self.mf.n_grow_vect())
imax = self._get_indices(min_box.big_end, 0)
return imax
def _fix_index(self, ii, imax, d):
"""Handle negative index, wrapping them as needed.
This depends on whether ghost cells are included. When true, the indices are
shifted by the number of ghost cells before being wrapped.
"""
nghosts = self._get_n_ghosts()
if self.include_ghosts:
ii += nghosts[d]
if ii < 0:
ii += imax
if self.include_ghosts:
ii -= nghosts[d]
return ii
def _find_start_stop(self, ii, imin, imax, d):
"""Given the input index, calculate the start and stop range of the indices.
Parameters
----------
ii: None, slice, integer
Input index, either None, a slice object, or an integer.
Note that ii can be negative.
imin: integer
The global lowest lower bound in the specified direction.
This can include the ghost cells.
imax: integer
The global highest upper bound in the specified direction.
This can include the ghost cells.
This should be the max index + 1.
d: integer
The dimension number, 0, 1, 2, or 3 (3 being the components)
If ii is a slice, the start and stop values are used directly,
unless they are None, then the lower or upper bound is used.
An assertion checks if the indices are within the bounds.
"""
if ii is None:
iistart = imin
iistop = imax
elif isinstance(ii, slice):
if ii.start is None:
iistart = imin
else:
iistart = self._fix_index(ii.start, imax, d)
if ii.stop is None:
iistop = imax
else:
iistop = self._fix_index(ii.stop, imax, d)
else:
ii = self._fix_index(ii, imax, d)
iistart = ii
iistop = ii + 1
assert imin <= iistart <= imax, Exception(f'Dimension {d+1} lower index is out of bounds')
assert imin <= iistop <= imax, Exception(f'Dimension {d+1} upper index is out of bounds')
return iistart, iistop
def _get_field(self, mfi):
"""Return the field at the given mfi.
If include ghosts is true, return the whole array, otherwise
return the interior slice that does not include the ghosts.
"""
# Note that the array will always have 4 dimensions.
# even when self.dim < 3.
# The transpose is taken since the Python array interface to Array4 in
# self.mf.array(mfi) is in C ordering.
# Note: transposing creates a view and not a copy.
device_arr4 = self.mf.array(mfi)
if cp is not None:
device_arr = cp.array(device_arr4, copy=False).T
else:
device_arr = np.array(device_arr4, copy=False).T
if not self.include_ghosts:
nghosts = self._get_n_ghosts()
device_arr = device_arr[tuple([slice(ng, -ng) for ng in nghosts[:self.dim]])]
return device_arr
def _get_intersect_slice(self, mfi, starts, stops, icstart, icstop):
"""Return the slices where the block intersects with the global slice.
If the block does not intersect, return None.
This also shifts the block slices by the number of ghost cells in the
MultiFab arrays since the arrays include the ghost cells.
Parameters
----------
mfi: MFIter
The MFIter instance for the current block,
starts: sequence
The minimum indices of the global slice.
These can be negative.
stops: sequence
The maximum indices of the global slice.
These can be negative.
icstart: integer
The minimum component index of the global slice.
These can be negative.
icstops: integer
The maximum component index of the global slice.
These can be negative.
Returns
-------
block_slices:
The slice of the intersection relative to the block
global_slices:
The slice of the intersection relative to the global array where the data from individual block will go
"""
box = mfi.tilebox()
if self.include_ghosts:
box.grow(self.mf.n_grow_vect())
ilo = self._get_indices(box.small_end, 0)
ihi = self._get_indices(box.big_end, 0)
# Add 1 to the upper end to be consistent with the slicing notation
ihi_p1 = [i + 1 for i in ihi]
i1 = np.maximum(starts, ilo)
i2 = np.minimum(stops, ihi_p1)
if np.all(i1 < i2):
block_slices = []
global_slices = []
for i in range(3):
block_slices.append(slice(i1[i] - ilo[i], i2[i] - ilo[i]))
global_slices.append(slice(i1[i] - starts[i], i2[i] - starts[i]))
block_slices.append(slice(icstart, icstop))
global_slices.append(slice(0, icstop - icstart))
return tuple(block_slices), tuple(global_slices)
else:
return None, None
def __getitem__(self, index):
"""Returns slice of the MultiFab using global indexing.
The shape of the object returned depends on the number of ix, iy and iz specified, which
can be from none to all three. Note that the values of ix, iy and iz are
relative to the fortran indexing, meaning that 0 is the lower boundary
of the whole domain, and in fortran ordering, i.e. [ix,iy,iz].
This allows negative indexing, though with ghosts cells included, the first n-ghost negative
indices will refer to the lower guard cells.
Parameters
----------
index: integer, or sequence of integers or slices, or Ellipsis
Index of the slice to return
"""
# Note that the index can have negative values (which wrap around) and has 1 added to the upper
# limit using python style slicing
if index == Ellipsis:
index = self.dim*[slice(None)]
elif isinstance(index, slice):
# If only one slice passed in, it was not wrapped in a list
index = [index]
if len(index) < self.dim+1:
# Add extra dims to index, including for the component.
# These are the dims left out and assumed to extend over the full size of the dim
index = list(index)
while len(index) < self.dim+1:
index.append(slice(None))
elif len(index) > self.dim+1:
raise Exception('Too many indices given')
# Expand the indices to length 3
ii = self._get_indices(index, None)
ic = index[-1]
# Global extent. These include the ghost cells when include_ghosts is True
ixmin, iymin, izmin = self._get_min_indices()
ixmax, iymax, izmax = self._get_max_indices()
# Setup the size of the array to be returned
ixstart, ixstop = self._find_start_stop(ii[0], ixmin, ixmax+1, 0)
iystart, iystop = self._find_start_stop(ii[1], iymin, iymax+1, 1)
izstart, izstop = self._find_start_stop(ii[2], izmin, izmax+1, 2)
icstart, icstop = self._find_start_stop(ic, 0, self.mf.n_comp(), 3)
# Gather the data to be included in a list to be sent to other processes
starts = [ixstart, iystart, izstart]
stops = [ixstop, iystop, izstop]
datalist = []
for mfi in self.mf:
block_slices, global_slices = self._get_intersect_slice(mfi, starts, stops, icstart, icstop)
if global_slices is not None:
# Note that the array will always have 4 dimensions.
device_arr = self._get_field(mfi)
if cp is not None:
# Copy the data from the device to the host
slice_arr = cp.asnumpy(device_arr[block_slices])
else:
slice_arr = device_arr[block_slices]
datalist.append((global_slices, slice_arr))
# Gather the data from all processors
if npes == 1:
all_datalist = [datalist]
else:
all_datalist = comm_world.allgather(datalist)
# Create the array to be returned
result_shape = (max(0, ixstop - ixstart),
max(0, iystop - iystart),
max(0, izstop - izstart),
max(0, icstop - icstart))
# Now, copy the data into the result array
result_global = None
for datalist in all_datalist:
for global_slices, f_arr in datalist:
if result_global is None:
# Delay allocation to here so that the type can be obtained
result_global = np.zeros(result_shape, dtype=f_arr.dtype)
result_global[global_slices] = f_arr
if result_global is None:
# Something went wrong with the index and no data was found. Return an empty array.
result_global = np.zeros(0)
# Remove dimensions of length 1, and if all dimensions
# are removed, return a scalar (that's what the [()] does)
return result_global.squeeze()[()]
def __setitem__(self, index, value):
"""Sets slices of a decomposed array.
The shape of the input object depends on the number of arguments specified, which can
be from none to all three.
This allows negative indexing, though with ghosts cells included, the first n-ghost negative
indices will refer to the lower guard cells.
Parameters
----------
index: integer, or sequence of integers or slices, or Ellipsis
The slice to set
value: scalar or array
Input value to assign to the specified slice of the MultiFab
"""
# Note that the index can have negative values (which wrap around) and has 1 added to the upper
# limit using python style slicing
if index == Ellipsis:
index = tuple(self.dim*[slice(None)])
elif isinstance(index, slice):
# If only one slice passed in, it was not wrapped in a list
index = [index]
if len(index) < self.dim+1:
# Add extra dims to index, including for the component.
# These are the dims left out and assumed to extend over the full size of the dim.
index = list(index)
while len(index) < self.dim+1:
index.append(slice(None))
elif len(index) > self.dim+1:
raise Exception('Too many indices given')
# Expand the indices to length 3
ii = self._get_indices(index, None)
ic = index[-1]
# Global extent. These include the ghost cells when include_ghosts is True
ixmin, iymin, izmin = self._get_min_indices()
ixmax, iymax, izmax = self._get_max_indices()
# Setup the size of the global array to be set
ixstart, ixstop = self._find_start_stop(ii[0], ixmin, ixmax+1, 0)
iystart, iystop = self._find_start_stop(ii[1], iymin, iymax+1, 1)
izstart, izstop = self._find_start_stop(ii[2], izmin, izmax+1, 2)
icstart, icstop = self._find_start_stop(ic, 0, self.mf.n_comp(), 3)
if isinstance(value, np.ndarray):
# Expand the shape of the input array to match the shape of the global array
# (it needs to be 4-D).
# This converts value to an array if needed, and the [...] grabs a view so
# that the shape change below doesn't affect value.
value3d = np.array(value)[...]
global_shape = list(value3d.shape)
# The shape of 1 is added for the extra dimensions and when index is an integer
# (in which case the dimension was not in the input array).
if not isinstance(ii[0], slice): global_shape[0:0] = [1]
if not isinstance(ii[1], slice): global_shape[1:1] = [1]
if not isinstance(ii[2], slice): global_shape[2:2] = [1]
if not isinstance(ic , slice) or len(global_shape) < 4: global_shape[3:3] = [1]
value3d.shape = global_shape
starts = [ixstart, iystart, izstart]
stops = [ixstop, iystop, izstop]
for mfi in self.mf:
block_slices, global_slices = self._get_intersect_slice(mfi, starts, stops, icstart, icstop)
if global_slices is not None:
mf_arr = self._get_field(mfi)
if isinstance(value, np.ndarray):
slice_value = value3d[global_slices]
if cp is not None:
# Copy data from host to device
slice_value = cp.asarray(value3d[global_slices])
mf_arr[block_slices] = slice_value
else:
mf_arr[block_slices] = value
def min(self, *args):
return self.mf.min(*args)
def max(self, *args):
return self.mf.max(*args)
def sum(self, *args):
return self.mf.sum(*args)
def min_index(self, *args):
return self.mf.minIndex(*args)
def max_index(self, *args):
return self.mf.maxIndex(*args)
def norm0(self, *args):
return self.mf.norm0(*args)
def ExWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Efield_aux[x]', level=level, include_ghosts=include_ghosts)
def EyWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Efield_aux[y]', level=level, include_ghosts=include_ghosts)
def EzWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Efield_aux[z]', level=level, include_ghosts=include_ghosts)
def BxWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Bfield_aux[x]', level=level, include_ghosts=include_ghosts)
def ByWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Bfield_aux[y]', level=level, include_ghosts=include_ghosts)
def BzWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Bfield_aux[z]', level=level, include_ghosts=include_ghosts)
def JxWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'current_fp[x]', level=level, include_ghosts=include_ghosts)
def JyWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'current_fp[y]', level=level, include_ghosts=include_ghosts)
def JzWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'current_fp[z]', level=level, include_ghosts=include_ghosts)
def ExFPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Efield_fp[x]', level=level, include_ghosts=include_ghosts)
def EyFPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Efield_fp[y]', level=level, include_ghosts=include_ghosts)
def EzFPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Efield_fp[z]', level=level, include_ghosts=include_ghosts)
def BxFPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Bfield_fp[x]', level=level, include_ghosts=include_ghosts)
def ByFPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Bfield_fp[y]', level=level, include_ghosts=include_ghosts)
def BzFPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Bfield_fp[z]', level=level, include_ghosts=include_ghosts)
def JxFPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'current_fp[x]', level=level, include_ghosts=include_ghosts)
def JyFPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'current_fp[y]', level=level, include_ghosts=include_ghosts)
def JzFPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'current_fp[z]', level=level, include_ghosts=include_ghosts)
def RhoFPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'rho_fp', level=level, include_ghosts=include_ghosts)
def PhiFPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'phi_fp', level=level, include_ghosts=include_ghosts)
def FFPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'F_fp', level=level, include_ghosts=include_ghosts)
def GFPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'G_fp', level=level, include_ghosts=include_ghosts)
def AxFPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'vector_potential_fp_nodal[x]', level=level, include_ghosts=include_ghosts)
def AyFPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'vector_potential_fp_nodal[y]', level=level, include_ghosts=include_ghosts)
def AzFPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'vector_potential_fp_nodal[z]', level=level, include_ghosts=include_ghosts)
def ExCPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Efield_cp[x]', level=level, include_ghosts=include_ghosts)
def EyCPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Efield_cp[y]', level=level, include_ghosts=include_ghosts)
def EzCPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Efield_cp[z]', level=level, include_ghosts=include_ghosts)
def BxCPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Bfield_cp[x]', level=level, include_ghosts=include_ghosts)
def ByCPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Bfield_cp[y]', level=level, include_ghosts=include_ghosts)
def BzCPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'Bfield_cp[z]', level=level, include_ghosts=include_ghosts)
def JxCPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'current_cp[x]', level=level, include_ghosts=include_ghosts)
def JyCPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'current_cp[y]', level=level, include_ghosts=include_ghosts)
def JzCPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'current_cp[z]', level=level, include_ghosts=include_ghosts)
def RhoCPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'rho_cp', level=level, include_ghosts=include_ghosts)
def FCPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'F_cp', level=level, include_ghosts=include_ghosts)
def GCPWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'G_cp', level=level, include_ghosts=include_ghosts)
def EdgeLengthsxWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'm_edge_lengths[x]', level=level, include_ghosts=include_ghosts)
def EdgeLengthsyWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'm_edge_lengths[y]', level=level, include_ghosts=include_ghosts)
def EdgeLengthszWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'm_edge_lengths[z]', level=level, include_ghosts=include_ghosts)
def FaceAreasxWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'm_face_areas[x]', level=level, include_ghosts=include_ghosts)
def FaceAreasyWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'm_face_areas[y]', level=level, include_ghosts=include_ghosts)
def FaceAreaszWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'm_face_areas[z]', level=level, include_ghosts=include_ghosts)
def JxFPAmpereWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'current_fp_ampere[x]', level=level, include_ghosts=include_ghosts)
def JyFPAmpereWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'current_fp_ampere[y]', level=level, include_ghosts=include_ghosts)
def JzFPAmpereWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name=f'current_fp_ampere[z]', level=level, include_ghosts=include_ghosts)
def ExFPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_E_fp[x]', level=level, include_ghosts=include_ghosts)
def EyFPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_E_fp[y]', level=level, include_ghosts=include_ghosts)
def EzFPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_E_fp[z]', level=level, include_ghosts=include_ghosts)
def BxFPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_B_fp[x]', level=level, include_ghosts=include_ghosts)
def ByFPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_B_fp[y]', level=level, include_ghosts=include_ghosts)
def BzFPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_B_fp[z]', level=level, include_ghosts=include_ghosts)
def JxFPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_j_fp[x]', level=level, include_ghosts=include_ghosts)
def JyFPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_j_fp[y]', level=level, include_ghosts=include_ghosts)
def JzFPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_j_fp[z]', level=level, include_ghosts=include_ghosts)
def FFPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_F_fp', level=level, include_ghosts=include_ghosts)
def GFPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_G_fp', level=level, include_ghosts=include_ghosts)
def ExCPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_E_cp[x]', level=level, include_ghosts=include_ghosts)
def EyCPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_E_cp[y]', level=level, include_ghosts=include_ghosts)
def EzCPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_E_cp[z]', level=level, include_ghosts=include_ghosts)
def BxCPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_B_cp[x]', level=level, include_ghosts=include_ghosts)
def ByCPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_B_cp[y]', level=level, include_ghosts=include_ghosts)
def BzCPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_B_cp[z]', level=level, include_ghosts=include_ghosts)
def JxCPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_j_cp[x]', level=level, include_ghosts=include_ghosts)
def JyCPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_j_cp[y]', level=level, include_ghosts=include_ghosts)
def JzCPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_j_cp[z]', level=level, include_ghosts=include_ghosts)
def FCPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_F_cp', level=level, include_ghosts=include_ghosts)
def GCPPMLWrapper(level=0, include_ghosts=False):
return _MultiFABWrapper(mf_name='pml_G_cp', level=level, include_ghosts=include_ghosts)
|