00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023 #ifndef BL_FABARRAY_H
00024 #define BL_FABARRAY_H
00025
00026
00027
00028
00029 #if defined(BL_OLD_STL)
00030 #include <string.h>
00031 #else
00032 #include <cstring>
00033 #endif
00034
00035 #include <limits>
00036 #include <map>
00037 #include <utility>
00038 #include <vector>
00039
00040 #include <BLassert.H>
00041 #include <PArray.H>
00042 #include <Array.H>
00043
00044 #include <Box.H>
00045 #include <BoxLib.H>
00046 #include <BoxArray.H>
00047 #include <BoxDomain.H>
00048 #include <FArrayBox.H>
00049 #include <DistributionMapping.H>
00050 #include <ParallelDescriptor.H>
00051 #include <ccse-mpi.H>
00052
00053
00054 class FabArrayBase
00055 {
00056 public:
00057
00058 FabArrayBase ();
00059 FabArrayBase (const BoxArray& bx, int nvar, int ngrow);
00060 FabArrayBase (const BoxArray& bx, int nvar, int ngrow, const DistributionMapping& map);
00061
00062 virtual ~FabArrayBase();
00063
00065
00066 int nGrow () const;
00067
00069
00070 int nComp () const;
00074 const BoxArray& boxArray () const;
00075
00079 virtual const Box& box (int K) const;
00080
00084 virtual Box fabbox (int K) const;
00085
00087
00088 int size () const;
00089
00091
00092 const DistributionMapping& DistributionMap () const;
00093
00094 protected:
00095
00096
00097
00098 mutable BoxArray boxarray;
00099 DistributionMapping distributionMap;
00100 int n_grow;
00101 int n_comp;
00102 };
00103
00104 class MFIter
00105 {
00106 public:
00107
00109
00110 explicit MFIter (const FabArrayBase& fabarray);
00111
00113
00114 const Box& validbox () const;
00115
00117
00118 Box fabbox () const;
00119
00121
00122 void operator++ ();
00123
00125
00126 bool isValid ();
00127
00129
00130 int index () const;
00131
00133
00134 const FabArrayBase& theFabArrayBase () const;
00135
00136 static void setDebugging (bool debugging);
00137
00138 protected:
00139
00140 static bool g_debugging;
00141
00142 const FabArrayBase& fabArray;
00143 int currentIndex;
00144 bool m_debugging;
00145 };
00146
00147
00148
00149
00150 template <class FAB> class FabArray;
00151 template <class FAB> class FabArrayCopyDescriptor;
00152
00153
00154
00156
00213 enum FabAlloc { Fab_noallocate = 0, Fab_allocate };
00214
00215 template <class FAB>
00216 class FabArray
00217 :
00218 public FabArrayBase
00219 {
00220 public:
00221
00222 typedef typename FAB::value_type value_type;
00223
00225
00226 FabArray ();
00227
00240 FabArray (const BoxArray& bxs,
00241 int nvar,
00242 int ngrow,
00243 FabAlloc mem_mode = Fab_allocate);
00244
00245 FabArray (const BoxArray& bxs,
00246 int nvar,
00247 int ngrow,
00248 const DistributionMapping& dm,
00249 FabAlloc mem_mode);
00250
00251
00253
00254 virtual ~FabArray ();
00255
00261 void define (const BoxArray& bxs,
00262 int nvar,
00263 int ngrow,
00264 FabAlloc mem_mode);
00265
00266 void define (const BoxArray& bxs,
00267 int nvar,
00268 int ngrow,
00269 const DistributionMapping& dm,
00270 FabAlloc mem_mode);
00271
00277 bool ok () const;
00278
00282 const FAB& operator[] (const MFIter& mfi) const;
00283
00284 const FAB& get (const MFIter& mfi) const;
00285
00287
00288 FAB& operator[] (const MFIter& mfi);
00289
00290 FAB& get (const MFIter& mfi);
00291
00295 const FAB& operator[] (int K) const;
00296
00298
00299 FAB& operator[] (int K);
00300
00302
00303 void setFab (int K, FAB* elem);
00304
00306
00307 void clear ();
00308
00310
00311 void setVal (value_type val);
00312 void operator= (const value_type& val);
00313
00318 void setVal (value_type val,
00319 int comp,
00320 int num_comp,
00321 int nghost = 0);
00322
00328 void setVal (value_type val,
00329 const Box& region,
00330 int comp,
00331 int num_comp,
00332 int nghost = 0);
00333
00337 void setVal (value_type val,
00338 int nghost);
00339
00344 void setVal (value_type val,
00345 const Box& region,
00346 int nghost);
00347
00349
00350 void setBndry (value_type val);
00351
00355 void setBndry (value_type val,
00356 int strt_comp,
00357 int ncomp);
00358
00364 void copy (const FabArray<FAB>& fa);
00365
00375 void copy (const FabArray<FAB>& src,
00376 int src_comp,
00377 int dest_comp,
00378 int num_comp);
00379
00383 void copy (FAB& dest) const;
00384
00389 void copy (FAB& dest,
00390 const Box& subbox) const;
00391
00397 void copy (FAB& dest,
00398 int src_comp,
00399 int dest_comp,
00400 int num_comp) const;
00401
00407 void copy (FAB& dest,
00408 const Box& subbox,
00409 int src_comp,
00410 int dest_comp,
00411 int num_comp) const;
00412
00415 void shift (const IntVect& v);
00416
00417 protected:
00418
00419 PArray<FAB> fabparray;
00420
00421 private:
00422
00423
00424
00425 FabArray (const FabArray<FAB>&);
00426 FabArray<FAB>& operator= (const FabArray<FAB>&);
00427
00428
00429
00430 void AllocFabs ();
00431 };
00432
00433
00434
00435
00436
00437 class FillBoxId
00438 {
00439 public:
00440
00441 FillBoxId ();
00442 FillBoxId (int newid, const Box& fillbox);
00443
00444 int Id () const;
00445 int FabIndex () const;
00446 void FabIndex (int fabindex);
00447 const Box& box () const;
00448
00449 private:
00450
00451 Box m_fillBox;
00452 int m_fillBoxId;
00453 int m_fabIndex;
00454 };
00455
00456
00457
00458
00459
00460 class CommDataCache
00461 {
00462 public:
00463
00464 CommDataCache ();
00465
00466 void operator= (const Array<CommData>& rhs);
00467
00468 bool isValid () const;
00469
00470 const Array<CommData>& theCommData () const;
00471
00472 private:
00473
00474 Array<CommData> m_commdata;
00475 bool m_valid;
00476 };
00477
00478 class FabArrayId
00479 {
00480 public:
00481
00482 explicit FabArrayId (int newid = -1);
00483
00484 int Id () const;
00485
00486 bool operator== (const FabArrayId& rhs) const;
00487
00488 private:
00489
00490 int fabArrayId;
00491 };
00492
00493
00494
00495
00496
00497
00498 enum FillType { FillLocally, FillRemotely, Unfillable };
00499
00500 template <class FAB>
00501 struct FabCopyDescriptor
00502 {
00503 FabCopyDescriptor ();
00504
00505 ~FabCopyDescriptor ();
00506
00507 FAB* localFabSource;
00508 Box subBox;
00509 int myProc;
00510 int copyFromProc;
00511 int copyFromIndex;
00512 int fillBoxId;
00513 int srcComp;
00514 int destComp;
00515 int nComp;
00516 FillType fillType;
00517 bool cacheDataAllocated;
00518
00519 private:
00520
00521
00522
00523 FabCopyDescriptor (const FabCopyDescriptor&);
00524 FabCopyDescriptor& operator= (const FabCopyDescriptor&);
00525 };
00526
00527 template <class FAB>
00528 FabCopyDescriptor<FAB>::FabCopyDescriptor ()
00529 :
00530 localFabSource(0),
00531 myProc(-1),
00532 copyFromProc(-1),
00533 copyFromIndex(-1),
00534 fillBoxId(-1),
00535 srcComp(-1),
00536 destComp(-1),
00537 nComp(-1),
00538 fillType(Unfillable),
00539 cacheDataAllocated(false)
00540 {}
00541
00542 template <class FAB>
00543 FabCopyDescriptor<FAB>::~FabCopyDescriptor ()
00544 {
00545 if (cacheDataAllocated)
00546 delete localFabSource;
00547 }
00548
00549
00550
00551
00552
00553
00554 template <class FAB>
00555 class FabArrayCopyDescriptor
00556 {
00557 typedef std::multimap<int,FabCopyDescriptor<FAB>*> FCDMap;
00558 typedef typename FCDMap::value_type FCDMapValueType;
00559 typedef typename FCDMap::iterator FCDMapIter;
00560 typedef typename FCDMap::const_iterator FCDMapConstIter;
00561
00562 public:
00563
00564 FabArrayCopyDescriptor ();
00565
00566 ~FabArrayCopyDescriptor ();
00567
00568 FabArrayId RegisterFabArray(FabArray<FAB> *fabarray);
00569
00570 FillBoxId AddBox (FabArrayId fabarrayid,
00571 const Box& destFabBox,
00572 BoxList* unfilledBoxes,
00573 bool bUseValidBox = true);
00574
00575 FillBoxId AddBox (FabArrayId fabarrayid,
00576 const Box& destFabBox,
00577 BoxList* unfilledBoxes,
00578 int srccomp,
00579 int destcomp,
00580 int numcomp,
00581 bool bUseValidBox = true);
00582
00583
00584
00585 FillBoxId AddBox (FabArrayId fabarrayid,
00586 const Box& destFabBox,
00587 BoxList* unfilledBoxes,
00588 int fabarrayindex,
00589 int srccomp,
00590 int destcomp,
00591 int numcomp,
00592 bool bUseValidBox = true);
00593
00594 void CollectData (Array<int>* snd_cache = 0,
00595 CommDataCache* cd_cache = 0);
00596
00597 void FillFab (FabArrayId fabarrayid,
00598 const FillBoxId& fillboxid,
00599 FAB& destFab);
00600
00601 void FillFab (FabArrayId fabarrayid,
00602 const FillBoxId& fillboxid,
00603 FAB& destFab,
00604 const Box& destBox);
00605
00606 void PrintStats () const;
00607
00608 bool DataAvailable () const { return dataAvailable; }
00609
00610 void clear ();
00611
00612 int nFabArrays () const { return fabArrays.size(); }
00613
00614 int nFabComTags () const { return fabComTagList.size(); }
00615
00616 int nFabCopyDescs () const { return fabCopyDescList.size(); }
00617
00618 protected:
00619
00620
00621
00622 void AddBoxDoIt (FabArrayId fabarrayid,
00623 const Box& destFabBox,
00624 BoxList* returnedUnfilledBoxes,
00625 int faindex,
00626 int srccomp,
00627 int destcomp,
00628 int numcomp,
00629 bool bUseValidBox,
00630 BoxDomain& unfilledBoxDomain,
00631 BoxList& filledBoxes);
00632
00633 std::vector<FabArray<FAB>*> fabArrays;
00634 std::vector<FCDMap> fabCopyDescList;
00635 std::vector<FabComTag> fabComTagList;
00636 int nextFillBoxId;
00637 bool dataAvailable;
00638
00639 private:
00640
00641
00642
00643 FabArrayCopyDescriptor (const FabArrayCopyDescriptor<FAB>&);
00644 FabArrayCopyDescriptor<FAB>& operator= (const FabArrayCopyDescriptor<FAB> &);
00645 };
00646
00647 inline
00648 int
00649 FabArrayBase::nGrow () const
00650 {
00651 return n_grow;
00652 }
00653
00654 inline
00655 const BoxArray&
00656 FabArrayBase::boxArray () const
00657 {
00658 return boxarray;
00659 }
00660
00661 inline
00662 const Box&
00663 FabArrayBase::box (int K) const
00664 {
00665 return boxarray[K];
00666 }
00667
00668 inline
00669 int
00670 FabArrayBase::size () const
00671 {
00672 return boxarray.size();
00673 }
00674
00675 inline
00676 int
00677 FabArrayBase::nComp () const
00678 {
00679 return n_comp;
00680 }
00681
00682 inline
00683 const DistributionMapping&
00684 FabArrayBase::DistributionMap () const
00685 {
00686 return distributionMap;
00687 }
00688
00689 inline
00690 const FabArrayBase&
00691 MFIter::theFabArrayBase () const
00692 {
00693 return fabArray;
00694 }
00695
00696 inline
00697 int
00698 MFIter::index () const
00699 {
00700 return currentIndex;
00701 }
00702
00703 inline
00704 int
00705 FillBoxId::Id () const
00706 {
00707 return m_fillBoxId;
00708 }
00709
00710 inline
00711 int
00712 FillBoxId::FabIndex () const
00713 {
00714 return m_fabIndex;
00715 }
00716
00717 inline
00718 void
00719 FillBoxId::FabIndex (int fabindex)
00720 {
00721 m_fabIndex = fabindex;
00722 }
00723
00724 inline
00725 const Box&
00726 FillBoxId::box () const
00727 {
00728 return m_fillBox;
00729 }
00730
00731 inline
00732 int
00733 FabArrayId::Id () const
00734 {
00735 return fabArrayId;
00736 }
00737
00738 template <class FAB>
00739 inline
00740 const FAB&
00741 FabArray<FAB>::operator[] (const MFIter& mfi) const
00742 {
00743 return fabparray[mfi.index()];
00744 }
00745
00746 template <class FAB>
00747 inline
00748 const FAB&
00749 FabArray<FAB>::get (const MFIter& mfi) const
00750 {
00751 return fabparray[mfi.index()];
00752 }
00753
00754 template <class FAB>
00755 inline
00756 FAB&
00757 FabArray<FAB>::operator[] (const MFIter& mfi)
00758 {
00759 return fabparray[mfi.index()];
00760 }
00761
00762 template <class FAB>
00763 inline
00764 FAB&
00765 FabArray<FAB>::get (const MFIter& mfi)
00766 {
00767 return fabparray[mfi.index()];
00768 }
00769
00770 template <class FAB>
00771 inline
00772 const FAB&
00773 FabArray<FAB>::operator[] (int K) const
00774 {
00775 return fabparray[K];
00776 }
00777
00778 template <class FAB>
00779 inline
00780 FAB&
00781 FabArray<FAB>::operator[] (int K)
00782 {
00783 return fabparray[K];
00784 }
00785
00786 template <class FAB>
00787 void
00788 FabArray<FAB>::clear ()
00789 {
00790 fabparray.clear();
00791 }
00792
00793 template <class FAB>
00794 void
00795 FabArray<FAB>::setVal (value_type val,
00796 int nghost)
00797 {
00798 setVal(val,0,n_comp,nghost);
00799 }
00800
00801 template <class FAB>
00802 void
00803 FabArray<FAB>::setVal (value_type val,
00804 const Box& region,
00805 int nghost)
00806 {
00807 setVal(val,region,0,n_comp,nghost);
00808 }
00809
00810 template <class FAB>
00811 FabArray<FAB>::FabArray ()
00812 :
00813 fabparray(0,PArrayManage)
00814 {}
00815
00816 template <class FAB>
00817 FabArray<FAB>::FabArray (const BoxArray& bxs,
00818 int nvar,
00819 int ngrow,
00820 FabAlloc alloc)
00821 :
00822 fabparray(0, PArrayManage)
00823 {
00824 define(bxs,nvar,ngrow,alloc);
00825 }
00826
00827 template <class FAB>
00828 FabArray<FAB>::FabArray (const BoxArray& bxs,
00829 int nvar,
00830 int ngrow,
00831 const DistributionMapping& dm,
00832 FabAlloc alloc)
00833 :
00834 fabparray(0, PArrayManage)
00835 {
00836 define(bxs,nvar,ngrow,dm,alloc);
00837 }
00838
00839
00840
00841
00842 template <class FAB>
00843 FabArray<FAB>::~FabArray ()
00844 {}
00845
00846 template <class FAB>
00847 bool
00848 FabArray<FAB>::ok () const
00849 {
00850 long isok = true;
00851
00852 for (MFIter fai(*this); fai.isValid() && isok; ++fai)
00853 {
00854 if (fabparray.defined(fai.index()))
00855 {
00856 if (get(fai).box() != BoxLib::grow(box(fai.index()),n_grow))
00857 {
00858 isok = false;
00859 }
00860 }
00861 else
00862 {
00863 isok = false;
00864 }
00865 }
00866
00867 ParallelDescriptor::ReduceLongAnd(isok);
00868
00869 return isok != 0;
00870 }
00871
00872 template <class FAB>
00873 void
00874 FabArray<FAB>::define (const BoxArray& bxs,
00875 int nvar,
00876 int ngrow,
00877 FabAlloc alloc)
00878 {
00879 BL_ASSERT(boxarray.size() == 0);
00880 n_grow = ngrow;
00881 n_comp = nvar;
00882 boxarray.define(bxs);
00883 distributionMap.define(boxarray,ParallelDescriptor::NProcsCFD());
00884 fabparray.resize(bxs.size());
00885 if (alloc == Fab_allocate)
00886 AllocFabs();
00887 }
00888
00889 template <class FAB>
00890 void
00891 FabArray<FAB>::define (const BoxArray& bxs,
00892 int nvar,
00893 int ngrow,
00894 const DistributionMapping& dm,
00895 FabAlloc alloc)
00896 {
00897 BL_ASSERT(boxarray.size() == 0);
00898 n_grow = ngrow;
00899 n_comp = nvar;
00900 boxarray.define(bxs);
00901 distributionMap = dm;
00902 fabparray.resize(bxs.size());
00903 if (alloc == Fab_allocate)
00904 AllocFabs();
00905 }
00906
00907 template <class FAB>
00908 void
00909 FabArray<FAB>::AllocFabs ()
00910 {
00911 for (MFIter fai(*this); fai.isValid(); ++fai)
00912 {
00913 Box tmp = BoxLib::grow(fai.validbox(), n_grow);
00914
00915 fabparray.set(fai.index(), new FAB(tmp, n_comp));
00916 }
00917 }
00918
00919 template <class FAB>
00920 void
00921 FabArray<FAB>::setFab (int boxno,
00922 FAB* elem)
00923 {
00924
00925
00926
00927 if (n_comp == 0)
00928 n_comp = elem->nComp();
00929
00930 BL_ASSERT(n_comp == elem->nComp());
00931 BL_ASSERT(boxarray.size() > 0);
00932 BL_ASSERT(elem->box() == BoxLib::grow(boxarray[boxno],n_grow));
00933 BL_ASSERT(!fabparray.defined(boxno));
00934 BL_ASSERT(distributionMap[boxno] == ParallelDescriptor::MyProc());
00935
00936 fabparray.set(boxno,elem);
00937 }
00938
00939 template <class FAB>
00940 void
00941 FabArray<FAB>::setBndry (value_type val)
00942 {
00943 setBndry(val, 0, n_comp);
00944 }
00945
00946 template <class FAB>
00947 void
00948 FabArray<FAB>::setBndry (value_type val,
00949 int strt_comp,
00950 int ncomp)
00951 {
00952 if (n_grow > 0)
00953 {
00954 for (MFIter fai(*this); fai.isValid(); ++fai)
00955 {
00956 get(fai).setComplement(val, fai.validbox(), strt_comp, ncomp);
00957 }
00958 }
00959 }
00960
00961 template <class FAB>
00962 void
00963 FabArray<FAB>::copy (const FabArray<FAB>& src,
00964 int scomp,
00965 int dcomp,
00966 int ncomp)
00967 {
00968
00969 if (boxarray == src.boxarray && distributionMap == src.distributionMap)
00970 {
00971 for (MFIter fai(*this); fai.isValid(); ++fai)
00972 {
00973 const Box& bx = fai.validbox();
00974 get(fai).copy(src[fai],bx,scomp,bx,dcomp,ncomp);
00975 }
00976
00977 return;
00978 }
00979
00980 const int MyProc = ParallelDescriptor::MyProc();
00981 const int NProcs = ParallelDescriptor::NProcs();
00982
00983 static Array<value_type*> fab_data(NProcs);
00984 static Array<int> indx(NProcs);
00985 static Array<MPI_Status> status(NProcs);
00986 static Array<MPI_Request> reqs(NProcs);
00987
00988 Array< std::vector<FabComTag> > SndTags(NProcs);
00989 Array< std::vector<FabComTag> > RcvTags(NProcs);
00990
00991 FabComTag tag;
00992
00993 for (int i = 0; i < size(); i++)
00994 {
00995 if (distributionMap[i] == MyProc)
00996 {
00997 for (int ii = 0; ii < src.boxarray.size(); ii++)
00998 {
00999 if (src.boxarray[ii].intersects(boxarray[i]))
01000 {
01001 Box bx = src.boxarray[ii] & boxarray[i];
01002
01003 if (src.distributionMap[ii] == MyProc)
01004 {
01005 fabparray[i].copy(src[ii],bx,scomp,bx,dcomp,ncomp);
01006 }
01007 else
01008 {
01009 tag.box = bx;
01010 tag.fabIndex = i;
01011
01012 RcvTags[src.distributionMap[ii]].push_back(tag);
01013 }
01014 }
01015 }
01016 }
01017 else
01018 {
01019 for (int ii = 0; ii < src.boxarray.size(); ii++)
01020 {
01021 if (src.distributionMap[ii] == MyProc)
01022 {
01023 if (src.boxarray[ii].intersects(boxarray[i]))
01024 {
01025 tag.box = src.boxarray[ii] & boxarray[i];
01026 tag.fabIndex = ii;
01027
01028 SndTags[distributionMap[i]].push_back(tag);
01029 }
01030 }
01031 }
01032 }
01033 }
01034
01035 if (NProcs == 1) return;
01036
01037 const int seqno = ParallelDescriptor::SeqNum();
01038
01039 int NWaits = 0;
01040
01041
01042
01043 for (int i = 0; i < NProcs; i++)
01044 {
01045 reqs[i] = MPI_REQUEST_NULL;
01046
01047 if (!RcvTags[i].empty())
01048 {
01049 NWaits++;
01050
01051 size_t N = 0;
01052
01053 for (unsigned int j = 0; j < RcvTags[i].size(); j++)
01054 N += RcvTags[i][j].box.numPts() * ncomp;
01055
01056 BL_ASSERT(N < std::numeric_limits<size_t>::max());
01057
01058 fab_data[i] = static_cast<value_type*>(BoxLib::The_Arena()->alloc(N*sizeof(value_type)));
01059
01060 reqs[i] = ParallelDescriptor::Arecv(fab_data[i],N,i,seqno).req();
01061 }
01062 }
01063
01064 FAB fab;
01065
01066
01067
01068 for (int i = 0; i < NProcs; i++)
01069 {
01070 if (!SndTags[i].empty())
01071 {
01072 size_t N = 0;
01073
01074 for (unsigned int j = 0; j < SndTags[i].size(); j++)
01075 N += SndTags[i][j].box.numPts() * ncomp;
01076
01077 BL_ASSERT(N < std::numeric_limits<size_t>::max());
01078
01079 value_type* data = static_cast<value_type*>(BoxLib::The_Arena()->alloc(N*sizeof(value_type)));
01080 value_type* dptr = data;
01081
01082 for (unsigned int j = 0; j < SndTags[i].size(); j++)
01083 {
01084 const Box& bx = SndTags[i][j].box;
01085 fab.resize(bx, ncomp);
01086 fab.copy(src[SndTags[i][j].fabIndex],bx,scomp,bx,0,ncomp);
01087 int count = bx.numPts() * ncomp;
01088 memcpy(dptr, fab.dataPtr(), count*sizeof(value_type));
01089 dptr += count;
01090 }
01091
01092 BL_ASSERT(data+N == dptr);
01093
01094 ParallelDescriptor::Send(data, N, i, seqno);
01095
01096 BoxLib::The_Arena()->free(data);
01097 }
01098 }
01099
01100
01101
01102 for (int completed; NWaits > 0; NWaits -= completed)
01103 {
01104 ParallelDescriptor::Waitsome(reqs, completed, indx, status);
01105
01106 for (int k = 0; k < completed; k++)
01107 {
01108 value_type* dptr = fab_data[indx[k]];
01109
01110 BL_ASSERT(!(dptr == 0));
01111
01112 for (unsigned int j = 0; j < RcvTags[indx[k]].size(); j++)
01113 {
01114 const Box& bx = RcvTags[indx[k]][j].box;
01115 fab.resize(bx, ncomp);
01116 int N = bx.numPts() * ncomp;
01117 memcpy(fab.dataPtr(), dptr, N*sizeof(value_type));
01118 fabparray[RcvTags[indx[k]][j].fabIndex].copy(fab,bx,0,bx,dcomp,ncomp);
01119 dptr += N;
01120 }
01121
01122 BoxLib::The_Arena()->free(fab_data[indx[k]]);
01123 }
01124 }
01125 }
01126
01127 template <class FAB>
01128 void
01129 FabArray<FAB>::copy (const FabArray<FAB>& src)
01130 {
01131 copy(src,0,0,nComp());
01132 }
01133
01134
01135
01136
01137
01138 template <class FAB>
01139 void
01140 FabArray<FAB>::copy (FAB& dest) const
01141 {
01142 copy(dest, dest.box(), 0, 0, dest.nComp());
01143 }
01144
01145 template <class FAB>
01146 void
01147 FabArray<FAB>::copy (FAB& dest,
01148 const Box& subbox) const
01149 {
01150 copy(dest, subbox, 0, 0, dest.nComp());
01151 }
01152
01153 template <class FAB>
01154 void
01155 FabArray<FAB>::copy (FAB& dest,
01156 int scomp,
01157 int dcomp,
01158 int ncomp) const
01159 {
01160 copy(dest, dest.box(), scomp, dcomp, ncomp);
01161 }
01162
01163 template <class FAB>
01164 void
01165 FabArray<FAB>::copy (FAB& dest,
01166 const Box& subbox,
01167 int scomp,
01168 int dcomp,
01169 int ncomp) const
01170 {
01171 BL_ASSERT(dcomp + ncomp <= dest.nComp());
01172
01173
01174 if (ParallelDescriptor::NProcs() == 1)
01175 {
01176 for (int j = 0; j < size(); ++j)
01177 {
01178 if (boxarray[j].intersects(subbox))
01179 {
01180 Box destbox = boxarray[j] & subbox;
01181
01182 dest.copy(fabparray[j],destbox,scomp,destbox,dcomp,ncomp);
01183 }
01184 }
01185
01186 return;
01187 }
01188
01189 FAB ovlp;
01190
01191 for (int i = 0; i < size(); i++)
01192 {
01193 if (subbox.intersects(boxarray[i]))
01194 {
01195 Box bx = subbox & boxarray[i];
01196
01197 ovlp.resize(bx,ncomp);
01198
01199 if (ParallelDescriptor::MyProc() == distributionMap[i])
01200 {
01201 ovlp.copy(fabparray[i],bx,scomp,bx,0,ncomp);
01202 }
01203
01204 const int N = bx.numPts()*ncomp;
01205
01206 ParallelDescriptor::Bcast(ovlp.dataPtr(),N,distributionMap[i]);
01207
01208 dest.copy(ovlp,bx,0,bx,dcomp,ncomp);
01209 }
01210 }
01211 }
01212
01213 template <class FAB>
01214 void
01215 FabArray<FAB>::setVal (value_type val)
01216 {
01217
01218 for (MFIter fai(*this); fai.isValid(); ++fai)
01219 {
01220 get(fai).setVal(val);
01221 }
01222 }
01223
01224 template <class FAB>
01225 inline
01226 void
01227 FabArray<FAB>::operator= (const value_type& val)
01228 {
01229 setVal(val);
01230 }
01231
01232 template <class FAB>
01233 void
01234 FabArray<FAB>::setVal (value_type val,
01235 int comp,
01236 int ncomp,
01237 int nghost)
01238 {
01239
01240 BL_ASSERT(nghost >= 0 && nghost <= n_grow);
01241 BL_ASSERT(comp+ncomp <= n_comp);
01242
01243 for (MFIter fai(*this); fai.isValid(); ++fai)
01244 {
01245 get(fai).setVal(val,BoxLib::grow(fai.validbox(),nghost), comp, ncomp);
01246 }
01247 }
01248
01249 template <class FAB>
01250 void
01251 FabArray<FAB>::setVal (value_type val,
01252 const Box& region,
01253 int comp,
01254 int ncomp,
01255 int nghost)
01256 {
01257
01258 BL_ASSERT(nghost >= 0 && nghost <= n_grow);
01259 BL_ASSERT(comp+ncomp <= n_comp);
01260
01261 for (MFIter fai(*this); fai.isValid(); ++fai)
01262 {
01263 Box b = BoxLib::grow(fai.validbox(),nghost) & region;
01264
01265 if (b.ok())
01266 get(fai).setVal(val, b, comp, ncomp);
01267 }
01268 }
01269
01270
01271 template <class FAB>
01272 void
01273 FabArray<FAB>::shift (const IntVect& v)
01274 {
01275
01276 for(int id(0); id < BL_SPACEDIM; ++id)
01277 {
01278 boxarray.shift(id, v[id]);
01279 }
01280 for (MFIter fai(*this); fai.isValid(); ++fai)
01281 {
01282 get(fai).shift(v);
01283 }
01284 }
01285
01286
01287 template <class FAB>
01288 FabArrayCopyDescriptor<FAB>::FabArrayCopyDescriptor ()
01289 :
01290 nextFillBoxId(0),
01291 dataAvailable(false)
01292 {}
01293
01294 template <class FAB>
01295 FabArrayId
01296 FabArrayCopyDescriptor<FAB>::RegisterFabArray(FabArray<FAB>* fabarray)
01297 {
01298 BL_ASSERT(fabArrays.size() == fabCopyDescList.size());
01299
01300 fabArrays.push_back(fabarray);
01301
01302 fabCopyDescList.resize(fabArrays.size(), FCDMap());
01303
01304 return FabArrayId(fabArrays.size() - 1);
01305 }
01306
01307 template <class FAB>
01308 void
01309 FabArrayCopyDescriptor<FAB>::AddBoxDoIt (FabArrayId fabarrayid,
01310 const Box& destFabBox,
01311 BoxList* returnedUnfilledBoxes,
01312 int faindex,
01313 int srccomp,
01314 int destcomp,
01315 int numcomp,
01316 bool bUseValidBox,
01317 BoxDomain& unfilledBoxDomain,
01318 BoxList& filledBoxes)
01319 {
01320 const int MyProc = ParallelDescriptor::MyProc();
01321
01322 FabArray<FAB>* fabArray = fabArrays[fabarrayid.Id()];
01323
01324 BL_ASSERT(faindex >= 0 && faindex < fabArray->size());
01325
01326 Box intersect = destFabBox;
01327
01328 if (bUseValidBox)
01329 {
01330 intersect &= fabArray->box(faindex);
01331 }
01332 else
01333 {
01334 intersect &= fabArray->fabbox(faindex);
01335 }
01336
01337 if (intersect.ok())
01338 {
01339 filledBoxes.push_back(intersect);
01340
01341 FabCopyDescriptor<FAB>* fcd = new FabCopyDescriptor<FAB>;
01342
01343 int remoteProc = fabArray->DistributionMap()[faindex];
01344 fcd->fillBoxId = nextFillBoxId;
01345 fcd->subBox = intersect;
01346 fcd->myProc = MyProc;
01347 fcd->copyFromProc = remoteProc;
01348 fcd->copyFromIndex = faindex;
01349 fcd->srcComp = srccomp;
01350 fcd->destComp = destcomp;
01351 fcd->nComp = numcomp;
01352
01353 if (MyProc == remoteProc)
01354 {
01355
01356
01357
01358 fcd->fillType = FillLocally;
01359 fcd->localFabSource = &(*fabArray)[faindex];
01360 }
01361 else
01362 {
01363
01364
01365
01366 FabComTag fabComTag;
01367
01368 dataAvailable = false;
01369 fcd->fillType = FillRemotely;
01370 fcd->localFabSource = new FAB(intersect, numcomp);
01371 fcd->cacheDataAllocated = true;
01372 fabComTag.fabArrayId = fabarrayid.Id();
01373 fabComTag.fillBoxId = nextFillBoxId;
01374 fabComTag.fabIndex = faindex;
01375 fabComTag.procThatNeedsData = MyProc;
01376 fabComTag.procThatHasData = remoteProc;
01377 fabComTag.box = intersect;
01378 fabComTag.srcComp = srccomp;
01379 fabComTag.destComp = destcomp;
01380 fabComTag.nComp = numcomp;
01381
01382
01383
01384 fabComTagList.push_back(fabComTag);
01385 }
01386
01387 fabCopyDescList[fabarrayid.Id()].insert(FCDMapValueType(fcd->fillBoxId,fcd));
01388
01389 if (!(returnedUnfilledBoxes == 0))
01390 {
01391 unfilledBoxDomain.rmBox(intersect);
01392 }
01393 }
01394 }
01395
01396 template <class FAB>
01397 FillBoxId
01398 FabArrayCopyDescriptor<FAB>::AddBox (FabArrayId fabarrayid,
01399 const Box& destFabBox,
01400 BoxList* returnedUnfilledBoxes,
01401 int srccomp,
01402 int destcomp,
01403 int numcomp,
01404 bool bUseValidBox)
01405 {
01406 BoxDomain unfilledBoxDomain(destFabBox.ixType());
01407 BoxList filledBoxes(destFabBox.ixType());
01408
01409 if (!(returnedUnfilledBoxes == 0))
01410 {
01411 unfilledBoxDomain.add(destFabBox);
01412 }
01413
01414 for (int i = 0, N = fabArrays[fabarrayid.Id()]->size(); i < N; i++)
01415 {
01416 AddBoxDoIt(fabarrayid,
01417 destFabBox,
01418 returnedUnfilledBoxes,
01419 i,
01420 srccomp,
01421 destcomp,
01422 numcomp,
01423 bUseValidBox,
01424 unfilledBoxDomain,
01425 filledBoxes);
01426 }
01427
01428 if (!(returnedUnfilledBoxes == 0))
01429 {
01430 returnedUnfilledBoxes->clear();
01431 (*returnedUnfilledBoxes) = unfilledBoxDomain.boxList();
01432 }
01433
01434 return FillBoxId(nextFillBoxId++, destFabBox);
01435 }
01436
01437 template <class FAB>
01438 FillBoxId
01439 FabArrayCopyDescriptor<FAB>::AddBox (FabArrayId fabarrayid,
01440 const Box& destFabBox,
01441 BoxList* returnedUnfilledBoxes,
01442 int fabarrayindex,
01443 int srccomp,
01444 int destcomp,
01445 int numcomp,
01446 bool bUseValidBox)
01447 {
01448 BoxDomain unfilledBoxDomain(destFabBox.ixType());
01449 BoxList filledBoxes(destFabBox.ixType());
01450
01451 if (!(returnedUnfilledBoxes == 0))
01452 {
01453 unfilledBoxDomain.add(destFabBox);
01454 }
01455
01456 AddBoxDoIt(fabarrayid,
01457 destFabBox,
01458 returnedUnfilledBoxes,
01459 fabarrayindex,
01460 srccomp,
01461 destcomp,
01462 numcomp,
01463 bUseValidBox,
01464 unfilledBoxDomain,
01465 filledBoxes);
01466
01467 if (!(returnedUnfilledBoxes == 0))
01468 {
01469 returnedUnfilledBoxes->clear();
01470 (*returnedUnfilledBoxes) = unfilledBoxDomain.boxList();
01471 }
01472
01473 return FillBoxId(nextFillBoxId++, destFabBox);
01474 }
01475
01476 template <class FAB>
01477 FillBoxId
01478 FabArrayCopyDescriptor<FAB>::AddBox (FabArrayId fabarrayid,
01479 const Box& destFabBox,
01480 BoxList* returnedUnfilledBoxes,
01481 bool bUseValidBox)
01482 {
01483 return AddBox(fabarrayid,
01484 destFabBox,
01485 returnedUnfilledBoxes,
01486 0,
01487 0,
01488 fabArrays[fabarrayid.Id()]->nComp(),
01489 bUseValidBox);
01490 }
01491
01492 template <class FAB>
01493 FabArrayCopyDescriptor<FAB>::~FabArrayCopyDescriptor()
01494 {
01495 clear();
01496 }
01497
01498 template <class FAB>
01499 void
01500 FabArrayCopyDescriptor<FAB>::clear ()
01501 {
01502 for (unsigned int i = 0; i < fabCopyDescList.size(); ++i)
01503 {
01504 FCDMapIter fmi = fabCopyDescList[i].begin();
01505
01506 for ( ; fmi != fabCopyDescList[i].end(); ++fmi)
01507 {
01508 delete (*fmi).second;
01509 }
01510 }
01511
01512 fabArrays.clear();
01513 fabCopyDescList.clear();
01514 fabComTagList.clear();
01515
01516 nextFillBoxId = 0;
01517 dataAvailable = false;
01518 }
01519
01520 template <class FAB>
01521 void
01522 FabArrayCopyDescriptor<FAB>::CollectData (Array<int>* snd_cache,
01523 CommDataCache* cd_cache)
01524 {
01525 typedef typename FAB::value_type value_type;
01526
01527 dataAvailable = true;
01528
01529 const int NProcs = ParallelDescriptor::NProcs();
01530
01531 if (NProcs == 1) return;
01532
01533
01534 const int MyProc = ParallelDescriptor::MyProc();
01535
01536 static Array<int> Snds(NProcs);
01537 static Array<int> Rcvs(NProcs);
01538 static Array<int> indx(NProcs);
01539 static Array<value_type*> fab_data(NProcs);
01540 static Array<MPI_Request> req_data(NProcs);
01541 static Array<MPI_Request> req_cd(NProcs);
01542 static Array<MPI_Status> status(NProcs);
01543
01544 Array<CommData> senddata;
01545 Array<CommData> recv_cd;
01546
01547 for (int i = 0; i < NProcs; i++)
01548 Snds[i] = Rcvs[i] = 0;
01549
01550 for (int i = 0; i < NProcs; i++)
01551 req_data[i] = req_cd[i] = MPI_REQUEST_NULL;
01552
01553 const int seqno_1 = ParallelDescriptor::SeqNum();
01554 const int seqno_2 = ParallelDescriptor::SeqNum();
01555
01556 int idx = 0, NumReqs = 0, NWaits;
01557
01558
01559
01560 for (unsigned int i = 0; i < fabComTagList.size(); i++)
01561 {
01562 BL_ASSERT(fabComTagList[i].box.ok());
01563 BL_ASSERT(fabComTagList[i].procThatNeedsData == MyProc);
01564 BL_ASSERT(fabComTagList[i].procThatHasData != MyProc);
01565
01566 Rcvs[fabComTagList[i].procThatHasData]++;
01567 }
01568 BL_ASSERT(Rcvs[MyProc] == 0);
01569
01570
01571
01572 #ifdef NDEBUG
01573 if (snd_cache == 0 || snd_cache->size() == 0)
01574 #endif
01575 {
01576 for (int i = 0; i < NProcs; i++)
01577 {
01578 ParallelDescriptor::Gather(&Rcvs[i], 1, Snds.dataPtr(), 1, i);
01579 }
01580
01581 BL_ASSERT(Snds[MyProc] == 0);
01582 }
01583
01584 if (snd_cache)
01585 {
01586 if (snd_cache->size() > 0)
01587 {
01588 BL_ASSERT(Snds == *snd_cache);
01589
01590 Snds = *snd_cache;
01591 }
01592 else
01593 {
01594 *snd_cache = Snds;
01595 }
01596 }
01597
01598 for (int i = 0; i < NProcs; i++)
01599 NumReqs += Snds[i];
01600
01601 recv_cd.resize(NumReqs);
01602
01603
01604
01605
01606 #ifdef NDEBUG
01607 if (cd_cache == 0 || !cd_cache->isValid())
01608 #endif
01609 {
01610
01611
01612
01613 BL_ASSERT(sizeof(CommData) == CommData::DIM*sizeof(int));
01614
01615 for (int i = 0; i < NProcs; i++)
01616 {
01617 if (Snds[i] > 0)
01618 {
01619 int* D = reinterpret_cast<int*>(&recv_cd[idx]);
01620 const int N = Snds[i] * CommData::DIM;
01621
01622 req_cd[i] = ParallelDescriptor::Arecv(D,N,i,seqno_1).req();
01623
01624 idx += Snds[i];
01625 }
01626 }
01627
01628 BL_ASSERT(idx == NumReqs);
01629
01630
01631
01632
01633 for (int k = 0, i = MyProc+1; k < NProcs; k++, i++)
01634 {
01635 i %= NProcs;
01636
01637 if (Rcvs[i] > 0)
01638 {
01639 senddata.resize(Rcvs[i]);
01640
01641 int Processed = 0;
01642
01643 for (unsigned int j = 0; j < fabComTagList.size(); j++)
01644 {
01645 if (fabComTagList[j].procThatHasData == i)
01646 {
01647 CommData data(0,
01648 fabComTagList[j].fabIndex,
01649 MyProc,
01650 0,
01651 fabComTagList[j].nComp,
01652 fabComTagList[j].srcComp,
01653 fabComTagList[j].fabArrayId,
01654 fabComTagList[j].box);
01655
01656 senddata[Processed++] = data;
01657 }
01658 }
01659
01660 BL_ASSERT(Processed == Rcvs[i]);
01661
01662 int* D = reinterpret_cast<int*>(senddata.dataPtr());
01663 const int N = senddata.size() * CommData::DIM;
01664
01665 ParallelDescriptor::Send(D, N, i, seqno_1);
01666 }
01667 }
01668
01669 NWaits = 0;
01670 for (int i = 0; i < NProcs; i++)
01671 if (req_cd[i] != MPI_REQUEST_NULL)
01672 NWaits++;
01673
01674 for (int completed; NWaits > 0; NWaits -= completed)
01675 {
01676 ParallelDescriptor::Waitsome(req_cd, completed, indx, status);
01677 }
01678 }
01679
01680 if (cd_cache)
01681 {
01682 if (cd_cache->isValid())
01683 {
01684 BL_ASSERT(recv_cd == cd_cache->theCommData());
01685
01686 recv_cd = cd_cache->theCommData();
01687 }
01688 else
01689 {
01690 *cd_cache = recv_cd;
01691 }
01692 }
01693
01694
01695
01696 for (int i = 0; i < NProcs; i++)
01697 {
01698 if (Rcvs[i] > 0)
01699 {
01700
01701
01702
01703 size_t N = 0;
01704
01705 for (unsigned int j = 0; j < fabComTagList.size(); j++)
01706 if (fabComTagList[j].procThatHasData == i)
01707 N += fabComTagList[j].box.numPts()*fabComTagList[j].nComp;
01708
01709 BL_ASSERT(N < std::numeric_limits<size_t>::max());
01710
01711 fab_data[i] = static_cast<value_type*>(BoxLib::The_Arena()->alloc(N*sizeof(value_type)));
01712
01713 req_data[i] = ParallelDescriptor::Arecv(fab_data[i],N,i,seqno_2).req();
01714 }
01715 }
01716
01717
01718
01719 idx = 0;
01720
01721 FAB fab;
01722
01723 for (int k = 0, i = MyProc+1; k < NProcs; k++, i++)
01724 {
01725 i %= NProcs;
01726
01727 int strt = 0;
01728 for (int j = 0; j < i; j++)
01729 strt += Snds[j];
01730
01731 if (Snds[i] > 0)
01732 {
01733 size_t N = 0;
01734
01735 for (int j = 0; j < Snds[i]; j++)
01736 N += recv_cd[strt+j].box().numPts() * recv_cd[strt+j].nComp();
01737
01738 BL_ASSERT(N < std::numeric_limits<size_t>::max());
01739
01740 value_type* data = static_cast<value_type*>(BoxLib::The_Arena()->alloc(N*sizeof(value_type)));
01741 value_type* dptr = data;
01742
01743 for (int j = 0; j < Snds[i]; j++)
01744 {
01745 const CommData& cd = recv_cd[strt+j];
01746
01747 BL_ASSERT(cd.id() == 0);
01748 BL_ASSERT(cd.fromproc() == i);
01749
01750 fab.resize(cd.box(), cd.nComp());
01751
01752 int count = fab.box().numPts() * fab.nComp();
01753
01754 fab.copy((*fabArrays[cd.fabarrayid()])[cd.fabindex()],
01755 fab.box(),
01756 cd.srcComp(),
01757 fab.box(),
01758 0,
01759 cd.nComp());
01760
01761 memcpy(dptr, fab.dataPtr(), count*sizeof(value_type));
01762
01763 dptr += count;
01764 }
01765
01766 BL_ASSERT(data+N == dptr);
01767
01768 ParallelDescriptor::Send(data, N, i, seqno_2);
01769
01770 BoxLib::The_Arena()->free(data);
01771
01772 idx += Snds[i];
01773 }
01774 }
01775
01776 BL_ASSERT(idx == NumReqs);
01777
01778
01779
01780 std::pair<FCDMapIter,FCDMapIter> match;
01781
01782 NWaits = 0;
01783 for (int i = 0; i < NProcs; i++)
01784 if (req_data[i] != MPI_REQUEST_NULL)
01785 NWaits++;
01786
01787 for (int completed; NWaits > 0; NWaits -= completed)
01788 {
01789 ParallelDescriptor::Waitsome(req_data, completed, indx, status);
01790
01791 for (int k = 0; k < completed; k++)
01792 {
01793 int Processed = 0;
01794 value_type* dptr = fab_data[indx[k]];
01795
01796 BL_ASSERT(!(dptr == 0));
01797
01798 for (unsigned int j = 0; j < fabComTagList.size(); j++)
01799 {
01800 if (fabComTagList[j].procThatHasData == indx[k])
01801 {
01802 const FabComTag& tag = fabComTagList[j];
01803
01804 match = fabCopyDescList[tag.fabArrayId].equal_range(tag.fillBoxId);
01805
01806 FCDMapIter fmi = match.first;
01807
01808 for ( ; fmi != match.second; ++fmi)
01809 {
01810 FabCopyDescriptor<FAB>* fcdp = (*fmi).second;
01811
01812 BL_ASSERT(fcdp->fillBoxId == tag.fillBoxId);
01813
01814 if (fcdp->subBox == tag.box)
01815 {
01816 value_type* dataPtr = fcdp->localFabSource->dataPtr();
01817 BL_ASSERT(!(dataPtr == 0));
01818 BL_ASSERT(fcdp->localFabSource->box() == tag.box);
01819 int N = tag.box.numPts()*tag.nComp;
01820 memcpy(dataPtr, dptr, N*sizeof(value_type));
01821 dptr += N;
01822 Processed++;
01823 break;
01824 }
01825 }
01826
01827 BL_ASSERT(!(fmi == match.second));
01828 }
01829 }
01830
01831 BL_ASSERT(Processed == Rcvs[indx[k]]);
01832
01833 BoxLib::The_Arena()->free(fab_data[indx[k]]);
01834
01835 fab_data[indx[k]] = 0;
01836 }
01837 }
01838 }
01839
01840 template <class FAB>
01841 void
01842 FabArrayCopyDescriptor<FAB>::FillFab (FabArrayId faid,
01843 const FillBoxId& fillboxid,
01844 FAB& destFab)
01845 {
01846 BL_ASSERT(dataAvailable);
01847
01848 std::pair<FCDMapIter,FCDMapIter> match = fabCopyDescList[faid.Id()].equal_range(fillboxid.Id());
01849
01850 for (FCDMapIter fmi = match.first; fmi != match.second; ++fmi)
01851 {
01852 FabCopyDescriptor<FAB>* fcdp = (*fmi).second;
01853
01854 BL_ASSERT(fcdp->fillBoxId == fillboxid.Id());
01855
01856 destFab.copy(*fcdp->localFabSource,
01857 fcdp->subBox,
01858 fcdp->fillType == FillLocally ? fcdp->srcComp : 0,
01859 fcdp->subBox,
01860 fcdp->destComp,
01861 fcdp->nComp);
01862 }
01863 }
01864
01865 template <class FAB>
01866 void
01867 FabArrayCopyDescriptor<FAB>::FillFab (FabArrayId faid,
01868 const FillBoxId& fillboxid,
01869 FAB& destFab,
01870 const Box& destBox)
01871 {
01872 BL_ASSERT(dataAvailable);
01873
01874 FCDMapIter fmi = fabCopyDescList[faid.Id()].lower_bound(fillboxid.Id());
01875
01876 BL_ASSERT(!(fmi == fabCopyDescList[faid.Id()].end()));
01877
01878 FabCopyDescriptor<FAB>* fcdp = (*fmi).second;
01879
01880 BL_ASSERT(fcdp->fillBoxId == fillboxid.Id());
01881
01882 BL_ASSERT(fcdp->subBox.sameSize(destBox));
01883
01884 destFab.copy(*fcdp->localFabSource,
01885 fcdp->subBox,
01886 fcdp->fillType == FillLocally ? fcdp->srcComp : 0,
01887 destBox,
01888 fcdp->destComp,
01889 fcdp->nComp);
01890
01891 BL_ASSERT(++fmi == fabCopyDescList[faid.Id()].upper_bound(fillboxid.Id()));
01892 }
01893
01894 template <class FAB>
01895 void
01896 FabArrayCopyDescriptor<FAB>::PrintStats () const
01897 {
01898 const int MyProc = ParallelDescriptor::MyProc();
01899
01900 std::cout << "----- "
01901 << MyProc
01902 << ": Parallel stats for FabArrayCopyDescriptor:" << '\n';
01903
01904 for (int fa = 0; fa < fabArrays.size(); ++fa)
01905 {
01906 std::cout << "fabArrays["
01907 << fa
01908 << "]->boxArray() = "
01909 << fabArrays[fa]->boxArray()
01910 << '\n';
01911 }
01912 }
01913
01914 #endif