Exemplo n.º 1
0
void
BoxList::catenate (BoxList& blist)
{
    BL_ASSERT(ixType() == blist.ixType());
    lbox.splice(lbox.end(), blist.lbox);
    BL_ASSERT(blist.isEmpty());
}
Exemplo n.º 2
0
const BoxList
BoxLib::intersect (const BoxList& bl,
                   const BoxList& br)
{
    BL_ASSERT(bl.ixType() == br.ixType());
    BoxList newbl(bl);
    return newbl.intersect(br);
}
Exemplo n.º 3
0
static
void
Boxes (const std::string& file,
       const std::string& name,
       BoxArray&          boxes,
       int&               boxesLevel)
{
    const std::string TheDflt = "default:";
    const std::string TheName = name + ":";

    std::ifstream is(file.c_str(),std::ios::in);

    if (!is.good())
        BoxLib::FileOpenFailed(file);

#define STRIP while( is.get() != '\n' )

    BoxArray ba_dflt;
    BoxArray ba_name;
    int bxLvl_dflt = -1;
    int bxLvl_name = -1;

    std::string line;

    while (std::getline(is,line))
    {
        if (line.empty() || line[0] == '#') continue;

        if (line == TheDflt || line == TheName)
        {
            bool dflt = (line == TheDflt) ? true : false;

            int N; int lvl; Box bx; BoxList bl;

            is >> N; STRIP;
            is >> lvl; STRIP;

            BL_ASSERT(N > 0);

            for (int i = 0; i < N; i++)
            {
                is >> bx; STRIP; bl.push_back(bx);
            }

			if (dflt)
			{
				ba_dflt = BoxArray(bl);
				bxLvl_dflt = lvl;
			}
			else
			{
				ba_name = BoxArray(bl);
				bxLvl_name = lvl;
			}
        }
    }
Exemplo n.º 4
0
bool
BoxList::operator== (const BoxList& rhs) const
{
    if ( !(size() == rhs.size()) ) return false;

    BoxList::const_iterator liter = begin(), riter = rhs.begin(), End = end();
    for (; liter != End; ++liter, ++riter)
        if ( !( *liter == *riter) )
            return false;
    return true;
}
Exemplo n.º 5
0
bool
BoxList::contains (const Box& b) const
{
    if (isEmpty()) return false;

    BL_ASSERT(ixType() == b.ixType());

    BoxList bnew = BoxLib::complementIn(b,*this);

    return bnew.isEmpty();
}
Exemplo n.º 6
0
void
BoxList::join (const BoxList& blist)
{
    BL_ASSERT(ixType() == blist.ixType());
    std::list<Box> lb = blist.lbox;
    lbox.splice(lbox.end(), lb);
}
Exemplo n.º 7
0
bool
BoxList::contains (const BoxList&  bl) const
{
    if (isEmpty() || bl.isEmpty()) return false;

    BL_ASSERT(ixType() == bl.ixType());

    if (!minimalBox().contains(bl.minimalBox())) return false;

    BoxArray ba(*this);

    for (const_iterator bli = bl.begin(), End = bl.end(); bli != End; ++bli)
        if (!ba.contains(*bli))
            return false;

    return true;
}
Exemplo n.º 8
0
const BoxList
BoxLib::complementIn (const Box&     b,
                      const BoxList& bl)
{
    BL_ASSERT(bl.ixType() == b.ixType());
    BoxList newb(b.ixType());
    newb.complementIn(b,bl);
    return newb;
}
Exemplo n.º 9
0
void Box::SaveBoxesToFile(const string& fileName, const BoxList& boxes)
{
    FILE* fout  = fopen(fileName.c_str(), "w");

    int n = boxes.size();
    fprintf(fout, "%d\n", n);
    for (int i = 0; i < n; i++)
        fprintf(fout, "%d %d %d\n", boxes[i].l, boxes[i].w, boxes[i].h);
    fclose(fout);
}
Exemplo n.º 10
0
BoxList Box::LoadBoxesFromFile(const string& fileName)
{
    FILE* fin  = fopen(fileName.c_str(), "r");
    BoxList boxes;
    int n;

    if (!fin)
    {
        printf("ERROR: NO such file '%s'\n", fileName.c_str());
        return boxes;
    }

    fscanf(fin, "%d", &n);
    for (int i = 0; i < n; i++)
    {
        int l, w, h;
        fscanf(fin, "%d%d%d", &l, &w, &h);
        boxes.push_back(Box(l, w, h));
    }
    fclose(fin);

    return boxes;
}
Exemplo n.º 11
0
void
AuxBoundaryData::initialize (const BoxArray& ba,
			     int             n_grow,
			     int             n_comp,
                             const Geometry& geom)
{
    BL_ASSERT(!m_initialized);

    const bool verbose   = false;
    const int  NProcs    = ParallelDescriptor::NProcs();
    const Real strt_time = ParallelDescriptor::second();

    m_ngrow = n_grow;

    BoxList gcells = BoxLib::GetBndryCells(ba,n_grow);
    //
    // Remove any intersections with periodically shifted valid region.
    //
    if (geom.isAnyPeriodic())
    {
        Box dmn = geom.Domain();

        for (int d = 0; d < BL_SPACEDIM; d++)
            if (!geom.isPeriodic(d)) 
                dmn.grow(d,n_grow);

        for (BoxList::iterator it = gcells.begin(); it != gcells.end(); )
        {
            const Box& isect = *it & dmn;

            if (isect.ok())
            {
                *it++ = isect;
            }
            else
            {
                gcells.remove(it++);
            }
        }
    }

    gcells.simplify();

    if (gcells.size() < NProcs)
    {
        gcells.maxSize(BL_SPACEDIM == 3 ? 64 : 128);
    }

    BoxArray nba(gcells);

    gcells.clear();

    if (nba.size() > 0)
    {
        m_fabs.define(nba, n_comp, 0, Fab_allocate);
    }
    else
    {
        m_empty = true;
    }

    if (verbose)
    {
        const int IOProc   = ParallelDescriptor::IOProcessorNumber();
        Real      run_time = ParallelDescriptor::second() - strt_time;
	const int sz       = nba.size();

#ifdef BL_LAZY
	Lazy::QueueReduction( [=] () mutable {
#endif
        ParallelDescriptor::ReduceRealMax(run_time,IOProc);
        if (ParallelDescriptor::IOProcessor()) 
            std::cout << "AuxBoundaryData::initialize() size = " << sz << ", time = " << run_time << '\n';
#ifdef BL_LAZY
	});
#endif
    }

    m_initialized = true;
}
Exemplo n.º 12
0
void ImagePane::paintEvent(QPaintEvent* /*event*/)
{
    QPainter painter(this);

    //draw the transformed version of the text-object image
    painter.drawImage(0, 0, transformedImage);

    //then we draw the bounding boxes
    //TODO inscriptions with at least one graph are a different color (?)

    //in case index numbers are visible, set font
    QFont font;
    font.setPixelSize(30/zoom); //TODO make font size dependent on resolution
    painter.setFont(font);
    QPen pen;
    pen.setWidth(0);
    pen.setColor(Qt::red);
    painter.setPen(pen);
    //make a list of bounding boxes according to current mode
    BoxList currentBoxList; //this is a list of all boxes
    currentBoxList.clear();
    switch(mode)
    {
    case SURFACE:
        currentBoxList.append(*surf); //list of one item, consting of the surface bounding box
        break;
    case INSCRIPTION:
        for(int i=0; i < surf->inscriptionCount(); i++)
            currentBoxList.insertBox(surf->inscrAt(i), i);
        break;
    case GRAPH:
        for(int i=0; i < surf->ptrInscrAt(currentInscrIndex)->count(); i++)
            currentBoxList.insertBox(surf->ptrInscrAt(currentInscrIndex)->at(i), i);
        break;
    default:
       break;
    }
    //iterate through the list of bounding boxes
    for (int i=0; i<currentBoxList.size(); i++)
    {
        BoundingBox currentBox = currentBoxList.at(i);

        //the bounding boxes need to be rotated and scaled
        QTransform boxTransform; //identity matrix
        //first we need to handle the bounding box's own rotation
        //by inverting the true matrix that was applied to the image
        //at the time each bounding box was created
        //(note: we don't need to worry about scale, as we took account of that
        //when the bounding box was created)
        boxTransform.rotate(currentBox.getRotation());
        boxTransform = QImage::trueMatrix(boxTransform,
                                          currentImage.width(), currentImage.height());
        boxTransform = boxTransform.inverted();
	
        //then we compound the above matrix with the current transformation of the image
        QTransform imageTrueTransform = QImage::trueMatrix(transform,
                                                           currentImage.width(), currentImage.height());
        painter.setWorldTransform(boxTransform * imageTrueTransform);
        //now draw the box
        //pen color is red; set the pen-color to green if this is the current box.
        if(i==currentBoxIndex)
        {
            QPen pen;
            pen.setWidth(0);
            pen.setColor(Qt::green);
            painter.setPen(pen);
        }
        painter.drawRect(currentBox);
        //and add an (optional) index number
        if(indexNumbersVisible)
        {
            painter.drawText(currentBox.left(), currentBox.top(), 50/zoom, 50/zoom,
                             Qt::AlignBottom,  QString("%1").arg(i+1)); //visible index, so base = 1, not zero
        }
        //return pen color to red (might be green)
        QPen pen;
        pen.setWidth(0);
        pen.setColor(Qt::red);
        painter.setPen(Qt::red);
    }

    //if label is not resized, it will stay large on zoom out
    //resulting in misleading / redundant scroll bars
    resize(transformedImage.size()); // keep label same size as image
}
Exemplo n.º 13
0
int
main (int   argc,
      char* argv[])
{
    BoxLib::Initialize(argc,argv);

    if (argc < 2)
        print_usage(argc,argv);

    ParmParse pp;

    if (pp.contains("help"))
        print_usage(argc,argv);

    bool verbose = false;
    pp.query("verbose",verbose);

    if (verbose>1)
        AmrData::SetVerbose(true);

    std::string infile; pp.get("infile",infile);
    std::string outfile_DEF;

    std::string outType = "tec";
#ifdef USE_TEC_BIN_IO
    bool doBin = true;
    pp.query("doBin",doBin);
    outfile_DEF = infile+(doBin ? ".plt" : ".dat" );
#else
    bool doBin=false;
    outfile_DEF = infile+".dat";
#endif

    // 
    bool connect_cc = true; pp.query("connect_cc",connect_cc);

    std::string outfile(outfile_DEF); pp.query("outfile",outfile);
    DataServices::SetBatchMode();
    Amrvis::FileType fileType(Amrvis::NEWPLT);
    DataServices dataServices(infile, fileType);

    if (!dataServices.AmrDataOk())
        DataServices::Dispatch(DataServices::ExitRequest, NULL);

    AmrData& amrData = dataServices.AmrDataRef();

    const Array<std::string>& names = amrData.PlotVarNames();

    Array<int> comps;
    if (int nc = pp.countval("comps"))
    {
        comps.resize(nc);
        pp.getarr("comps",comps,0,nc);
    }
    else
    {
        int sComp = 0;
        pp.query("sComp",sComp);
        int nComp = amrData.NComp();
        pp.query("nComp",nComp);
        BL_ASSERT(sComp+nComp <= amrData.NComp());
        comps.resize(nComp);
        for (int i=0; i<nComp; ++i)
            comps[i] = sComp + i;
    }

    Box subbox;
    if (int nx=pp.countval("box"))
    {
        Array<int> barr;
        pp.getarr("box",barr,0,nx);
        int d=BL_SPACEDIM;
        BL_ASSERT(barr.size()==2*d);
        subbox=Box(IntVect(D_DECL(barr[0],barr[1],barr[2])),
                   IntVect(D_DECL(barr[d],barr[d+1],barr[d+2]))) & amrData.ProbDomain()[0];

    } else {

        subbox = amrData.ProbDomain()[0];
    }

    int finestLevel = amrData.FinestLevel(); pp.query("finestLevel",finestLevel);
    int Nlev = finestLevel + 1;
    Array<BoxArray> gridArray(Nlev);
    Array<Box> subboxArray(Nlev);

    int nGrowPer = 0; pp.query("nGrowPer",nGrowPer);
    PArray<Geometry> geom(Nlev);

    Array<Real> LO(BL_SPACEDIM,0);
    Array<Real> HI(BL_SPACEDIM,1);
    RealBox rb(LO.dataPtr(),HI.dataPtr());
    int coordSys = 0;
    Array<int> isPer(BL_SPACEDIM,0);

    for (int lev=0; lev<Nlev; ++lev)
    {
        subboxArray[lev]
            = (lev==0 ? subbox
               : Box(subboxArray[lev-1]).refine(amrData.RefRatio()[lev-1]));

        geom.set(lev,new Geometry(amrData.ProbDomain()[lev],&rb,coordSys,const_cast<int*>(isPer.dataPtr())));

        if (nGrowPer>0 && lev==0)
        {
            for (int i=0; i<BL_SPACEDIM; ++i)
            {
                if (geom[lev].isPeriodic(i))
                {
                    if (subboxArray[lev].smallEnd()[i] == amrData.ProbDomain()[lev].smallEnd()[i])
                        subboxArray[lev].growLo(i,nGrowPer);
                    if (subboxArray[lev].bigEnd()[i] == amrData.ProbDomain()[lev].bigEnd()[i])
                        subboxArray[lev].growHi(i,nGrowPer);
                }
            }
        }

        gridArray[lev] = BoxLib::intersect(amrData.boxArray(lev), subboxArray[lev]);

        if (nGrowPer>0 && geom[lev].isAnyPeriodic() && gridArray[lev].size()>0)
        {
	  //const Box& probDomain = amrData.ProbDomain()[lev];
            const BoxArray& ba = amrData.boxArray(lev);
            BoxList bladd;
            Array<IntVect> shifts;
            for (int i=0; i<ba.size(); ++i)
            {
                geom[lev].periodicShift(subboxArray[lev],ba[i],shifts);
                for (int j=0; j<shifts.size(); ++j)
                {
                    Box ovlp = Box(ba[i]).shift(shifts[j]) & subboxArray[lev];
                    if (ovlp.ok())
                        bladd.push_back(ovlp);
                }
            }
            bladd.simplify();
            BoxList blnew(gridArray[lev]);
            blnew.join(bladd);
            gridArray[lev] = BoxArray(blnew);
        }

        if (!gridArray[lev].size())
        {
            Nlev = lev;
            finestLevel = Nlev-1;
            gridArray.resize(Nlev);
            subboxArray.resize(Nlev);
        }
    }

    const int nGrow = 1;
    typedef BaseFab<Node> NodeFab;
    typedef FabArray<NodeFab> MultiNodeFab; 
    PArray<MultiNodeFab> nodes(Nlev,PArrayManage);

    std::cerr << "Before nodes allocated" << endl;
    for (int lev=0; lev<Nlev; ++lev)
        nodes.set(lev,new MultiNodeFab(gridArray[lev],1,nGrow));
    std::cerr << "After nodes allocated" << endl;

    int cnt = 0;
    typedef std::map<Node,int> NodeMap;
    NodeMap nodeMap;
    for (int lev=0; lev<Nlev; ++lev)
    {
        for (MFIter fai(nodes[lev]); fai.isValid(); ++fai)
        {
            NodeFab& ifab = nodes[lev][fai];
            const Box& box = ifab.box() & subboxArray[lev];
            for (IntVect iv=box.smallEnd(); iv<=box.bigEnd(); box.next(iv))
                ifab(iv,0) = Node(iv,lev,fai.index(),Node::VALID);
        }
            
        if (lev != 0)
        {
            const int ref = amrData.RefRatio()[lev-1];
            const Box& rangeBox = Box(IntVect::TheZeroVector(),
                                     (ref-1)*IntVect::TheUnitVector());

            BoxArray bndryCells = GetBndryCells(nodes[lev].boxArray(),ref,geom[lev]);

            for (MFIter fai(nodes[lev]); fai.isValid(); ++fai)
            {
                const Box& box = BoxLib::grow(fai.validbox(),ref) & subboxArray[lev];
                NodeFab& ifab = nodes[lev][fai];
                std::vector< std::pair<int,Box> > isects = bndryCells.intersections(box);
                for (int i = 0; i < isects.size(); i++)
                {
                    Box co = isects[i].second & fai.validbox();
                    if (co.ok())
                        std::cout << "BAD ISECTS: " << co << std::endl;

                    const Box& dstBox = isects[i].second;
                    const Box& srcBox = BoxLib::coarsen(dstBox,ref);

                    NodeFab dst(dstBox,1);
                    for (IntVect iv(srcBox.smallEnd());
                         iv<=srcBox.bigEnd();
                         srcBox.next(iv))
                    {
                        const IntVect& baseIV = ref*iv;
                        for (IntVect ivt(rangeBox.smallEnd());ivt<=rangeBox.bigEnd();rangeBox.next(ivt))
                            dst(baseIV + ivt,0) = Node(iv,lev-1,-1,Node::VALID);
                    }
                    const Box& ovlp = dstBox & ifab.box();

                    Box mo = ovlp & fai.validbox();
                    if (mo.ok())
                    {
                        std::cout << "BAD OVERLAP: " << mo << std::endl;
                        std::cout << "         vb: " << fai.validbox() << std::endl;
                    }
                    if (ovlp.ok())
                        ifab.copy(dst,ovlp,0,ovlp,0,1);
                }
            }
        }

        // Block out cells covered by finer grid
        if (lev < finestLevel)
        {
            const BoxArray coarsenedFineBoxes =
                BoxArray(gridArray[lev+1]).coarsen(amrData.RefRatio()[lev]);
                
            for (MFIter fai(nodes[lev]); fai.isValid(); ++fai)
            {
                NodeFab& ifab = nodes[lev][fai];
                const Box& box = ifab.box();
                std::vector< std::pair<int,Box> > isects = coarsenedFineBoxes.intersections(box);
                for (int i = 0; i < isects.size(); i++)
                {
                    const Box& ovlp = isects[i].second;
                    for (IntVect iv=ovlp.smallEnd(); iv<=ovlp.bigEnd(); ovlp.next(iv))
                        ifab(iv,0) = Node(iv,lev,fai.index(),Node::COVERED);
                }
            }
        }

        // Add the unique nodes from this level to the list
        for (MFIter fai(nodes[lev]); fai.isValid(); ++fai)
        {
            NodeFab& ifab = nodes[lev][fai];
            const Box& box = fai.validbox() & subboxArray[lev];
            for (IntVect iv(box.smallEnd()); iv<=box.bigEnd(); box.next(iv))
            {
                if (ifab(iv,0).type == Node::VALID)
                {
                    if (ifab(iv,0).level != lev) 
                        std::cout << "bad level: " << ifab(iv,0) << std::endl;
                    nodeMap[ifab(iv,0)] = cnt++;
                }
            }
        }
    }

    std::cerr << "After nodeMap built, size=" << nodeMap.size() << endl;

    typedef std::set<Element> EltSet;
    EltSet elements;
    for (int lev=0; lev<Nlev; ++lev)
    {
        for (MFIter fai(nodes[lev]); fai.isValid(); ++fai)
        {
            NodeFab& ifab = nodes[lev][fai];                        
            Box box = ifab.box() & subboxArray[lev];

            for (int dir=0; dir<BL_SPACEDIM; ++dir)
                box.growHi(dir,-1);

            for (IntVect iv(box.smallEnd()); iv<=box.bigEnd(); box.next(iv))
            {
#if (BL_SPACEDIM == 2)
                const Node& n1 = ifab(iv,0);
                const Node& n2 = ifab(IntVect(iv).shift(BoxLib::BASISV(0)),0);
                const Node& n3 = ifab(IntVect(iv).shift(IntVect::TheUnitVector()),0);
                const Node& n4 = ifab(IntVect(iv).shift(BoxLib::BASISV(1)),0);
                if (n1.type==Node::VALID && n2.type==Node::VALID &&
                    n3.type==Node::VALID && n4.type==Node::VALID )
                    elements.insert(Element(n1,n2,n3,n4));
#else
                const IntVect& ivu = IntVect(iv).shift(BoxLib::BASISV(2));
                const Node& n1 = ifab(iv ,0);
                const Node& n2 = ifab(IntVect(iv ).shift(BoxLib::BASISV(0)),0);
                const Node& n3 = ifab(IntVect(iv ).shift(BoxLib::BASISV(0)).shift(BoxLib::BASISV(1)),0);
                const Node& n4 = ifab(IntVect(iv ).shift(BoxLib::BASISV(1)),0);
                const Node& n5 = ifab(ivu,0);
                const Node& n6 = ifab(IntVect(ivu).shift(BoxLib::BASISV(0)),0);
                const Node& n7 = ifab(IntVect(ivu).shift(BoxLib::BASISV(0)).shift(BoxLib::BASISV(1)),0);
                const Node& n8 = ifab(IntVect(ivu).shift(BoxLib::BASISV(1)),0);
                if (n1.type==Node::VALID && n2.type==Node::VALID &&
                    n3.type==Node::VALID && n4.type==Node::VALID &&
                    n5.type==Node::VALID && n6.type==Node::VALID &&
                    n7.type==Node::VALID && n8.type==Node::VALID )
                    elements.insert(Element(n1,n2,n3,n4,n5,n6,n7,n8));
#endif
            }
        }
    }

    int nElts = (connect_cc ? elements.size() : nodeMap.size());
    std::cerr << "Before connData allocated " << elements.size() << " elements"  << endl;
    Array<int> connData(MYLEN*nElts);
    std::cerr << "After connData allocated " << elements.size() << " elements" << endl;

    if (connect_cc) {
        cnt = 0;
        for (EltSet::const_iterator it = elements.begin(); it!=elements.end(); ++it)
        {
            for (int j=0; j<MYLEN; ++j)
            {
                const NodeMap::const_iterator noit = nodeMap.find(*((*it).n[j]));
                if (noit == nodeMap.end())
                {
                    std::cout << "Node not found in node map" << std::endl;
                    std::cout << *((*it).n[j]) << std::endl;
                }
                else
                {
                    connData[cnt++] = noit->second+1;
                }
            }
        }
    } else {
        cnt = 1;
        for (int i=0; i<nElts; ++i) {
            for (int j=0; j<MYLEN; ++j) {
                connData[i*MYLEN+j] = cnt++;
            }
        }
    }

    std::cerr << "Final elements built" << endl;

    // Invert the map
    std::vector<Node> nodeVect(nodeMap.size());
    for (NodeMap::const_iterator it=nodeMap.begin(); it!=nodeMap.end(); ++it)
    {
        if (it->second>=nodeVect.size() || it->second<0)
            std::cout << "Bad id: " << it->second << "  bad node: " << it->first << std::endl;
        BL_ASSERT(it->second>=0);
        BL_ASSERT(it->second<nodeVect.size());
        nodeVect[it->second] = (*it).first;
    }
    std::cerr << "Final nodeVect built (" << nodeVect.size() << " nodes)" << endl;
		
    nodeMap.clear();
    elements.clear();
    nodes.clear();

    std::cerr << "Temp nodes, elements cleared" << endl;

    PArray<MultiFab> fileData(Nlev);
    int ng = nGrowPer;
    for (int lev=0; lev<Nlev; ++lev)
    {
        if (lev!=0)
            ng *= amrData.RefRatio()[lev-1];
        const BoxArray& ba = gridArray[lev];
        fileData.set(lev,new MultiFab(ba,comps.size(),0));
        fileData[lev].setVal(1.e30);
        std::cerr << "My data set alloc'd at lev=" << lev << endl;

        MultiFab pData, pDataNG;
        if (ng>0 && geom[lev].isAnyPeriodic())
        {
            const Box& pd = amrData.ProbDomain()[lev];
            //const BoxArray& vba = amrData.boxArray(lev);
            Box shrunkenDomain = pd;
            for (int i=0; i<BL_SPACEDIM; ++i)
                if (geom[lev].isPeriodic(i))
                    shrunkenDomain.grow(i,-ng);

            const BoxArray edgeVBoxes = BoxLib::boxComplement(pd,shrunkenDomain);
            pData.define(edgeVBoxes,1,ng,Fab_allocate);
            pDataNG.define(BoxArray(edgeVBoxes).grow(ng),1,0,Fab_allocate);
        }

        for (int i=0; i<comps.size(); ++i)
        {
            BoxArray tmpBA = BoxLib::intersect(fileData[lev].boxArray(),amrData.ProbDomain()[lev]);
            MultiFab tmpMF(tmpBA,1,0);
            tmpMF.setVal(2.e30);
            amrData.FillVar(tmpMF,lev,names[comps[i]],0);
            fileData[lev].copy(tmpMF,0,i,1);
            if (ng>0 && geom[lev].isAnyPeriodic())
            {
                pData.setVal(3.e30);
                pDataNG.copy(tmpMF);
                for (MFIter mfi(pData); mfi.isValid(); ++mfi)
                    pData[mfi].copy(pDataNG[mfi]);
                amrData.FillVar(pData,lev,names[comps[i]],0);
                geom[lev].FillPeriodicBoundary(pData);
                for (MFIter mfi(pData); mfi.isValid(); ++mfi)
                    pDataNG[mfi].copy(pData[mfi]);
                fileData[lev].copy(pDataNG,0,i,1);
            }            
        }

        if (fileData[lev].max(0) > 1.e29)
        {
            std::cerr << "Bad mf data" << std::endl;
            VisMF::Write(fileData[lev],"out.mfab");
            BoxLib::Abort();
        }
    }
    std::cerr << "File data loaded" << endl;

    int nNodesFINAL = (connect_cc  ?  nodeVect.size() : nElts*MYLEN );
    int nCompsFINAL = BL_SPACEDIM+comps.size();
    FABdata tmpData(nNodesFINAL,nCompsFINAL);
    int tmpDatLen = nCompsFINAL*nNodesFINAL;
    std::cerr << "Final node data allocated (size=" << tmpDatLen << ")" << endl;

    //const Array<Array<Real> >& dxLevel = amrData.DxLevel();  // Do not trust dx from file...compute our own
    Array<Array<Real> > dxLevel(Nlev);
    for (int i=0; i<Nlev; ++i)
    {
        dxLevel[i].resize(BL_SPACEDIM);
        for (int j=0; j<BL_SPACEDIM; ++j)
            dxLevel[i][j] = amrData.ProbSize()[j]/amrData.ProbDomain()[i].length(j);
    }
    const Array<Real>& plo = amrData.ProbLo();

    // Handy structures for loading data in usual fab ordering (transpose of mef/flt ordering)
#define BIN_POINT
#undef BIN_POINT
#ifdef BIN_POINT
    Real* data = tmpData.dataPtr();
#else /* BLOCK ordering */
    Array<Real*> fdat(comps.size()+BL_SPACEDIM);
    for (int i=0; i<fdat.size(); ++i)
        fdat[i] = tmpData.dataPtr(i);
#endif

    cnt = 0;
    int Nnodes = nodeVect.size();
    int jGridPrev = 0;
    int levPrev = -1;
    for (int i=0; i<Nnodes; ++i)
    {
        const Node& node = nodeVect[i];

        const Array<Real>& dx = dxLevel[node.level];
        const IntVect& iv = node.iv;

        const BoxArray& grids = fileData[node.level].boxArray();
        int jGrid = node.grid;
        if (jGrid<0)
        {
            bool found_it = false;
            
            // Try same grid as last time, otherwise search list
            if (node.level==levPrev && grids[jGridPrev].contains(iv))
            {
                found_it = true;
                jGrid = jGridPrev;
            }
            for (int j=0; j<grids.size() && !found_it; ++j)
            {
                if (grids[j].contains(iv))
                {
                    found_it = true;
                    jGrid = j;
                }
            }
            BL_ASSERT(found_it);
        }
	// Remember these for next time
	levPrev = node.level;
	jGridPrev = jGrid;

	Array<IntVect> ivt;
        if (connect_cc) {
            ivt.resize(1,iv);
        } else {
            ivt.resize(D_PICK(1,4,8),iv);
            ivt[1] += BoxLib::BASISV(0);
            ivt[2] = ivt[1] + BoxLib::BASISV(1);
            ivt[3] += BoxLib::BASISV(1);
#if BLSPACEDIM==3
            for (int n=0; n<4; ++n) {
                ivt[4+n] = iv[n] + BoxLib::BASISV(2);
            }
#endif
        }

	for (int j=0; j<ivt.size(); ++j) {

            Real offset = (connect_cc ? 0.5 : 0);

#ifdef BIN_POINT
	  for (int dir=0; dir<BL_SPACEDIM; ++dir)
            data[cnt++] = plo[dir] + (ivt[j][dir] + offset) * dx[dir];
#else /* BLOCK ordering */
	  for (int dir=0; dir<BL_SPACEDIM; ++dir)
            fdat[dir][cnt] = plo[dir] + (ivt[j][dir] + offset) * dx[dir];
#endif /* BIN_POINT */


#ifdef BIN_POINT
	  for (int n=0; n<comps.size(); ++n) {
            data[cnt++] = fileData[node.level][jGrid](iv,n);
          }
#else /* BLOCK ordering */
	  for (int n=0; n<comps.size(); ++n) {
              fdat[n+BL_SPACEDIM][cnt] = fileData[node.level][jGrid](iv,n);
          }
	  cnt++;
#endif /* BIN_POINT */
	}
    }

    //Collate(tmpData,connData,MYLEN);

    //
    // Write output
    //
    const int nState = BL_SPACEDIM + comps.size();
    std::string vars = D_TERM("X"," Y"," Z");
    for (int j=0; j<comps.size(); ++j)
        vars += " " + amrData.PlotVarNames()[comps[j]];

    if (outType == "tec")
    {

#ifdef BIN_POINT
        string block_or_point = "FEPOINT";
#else /* BLOCK ordering */
        string block_or_point = "FEBLOCK";
#endif


        if (doBin)
        {
#ifdef USE_TEC_BIN_IO
            INTEGER4 Debug = 0;
            INTEGER4 VIsDouble = 1;
            INTEGER4 EltID = D_PICK(0,1,3);
            TECINI((char*)"Pltfile data", (char*)vars.c_str(), (char*)outfile.c_str(), (char*)".", &Debug, &VIsDouble);
            INTEGER4 nPts = nNodesFINAL;
            TECZNE((char*)infile.c_str(), &nPts, &nElts, &EltID, (char*)block_or_point.c_str(), NULL);
            TECDAT(&tmpDatLen,tmpData.fab.dataPtr(),&VIsDouble);
            TECNOD(connData.dataPtr());
            TECEND();
#else
            BoxLib::Abort("Need to recompile with USE_TEC_BIN_IO defined");
#endif
        }
        else
        {
            std::ofstream osf(outfile.c_str(),std::ios::out);
            osf << D_TERM("VARIABLES= \"X\"", " \"Y\"", " \"Z\"");
            for (int j=0; j<comps.size(); ++j)
                osf << " \""  << amrData.PlotVarNames()[comps[j]] << "\"";
            char buf[100];
            sprintf(buf,"%g",amrData.Time());
            osf << endl << "ZONE T=\"" << infile << " time = " << buf
                << "\", N=" << nNodesFINAL << ", E=" << nElts << ", F=" << "FEPOINT" << " ET="
                //<< "\", N=" << nPts << ", E=" << nElts << ", F=" << block_or_point << " ET="
                << D_PICK("POINT","QUADRILATERAL","BRICK") << endl;
            
            for (int i=0; i<nNodesFINAL; ++i)
            {
                for (int j=0; j<nState; ++j)
                    osf << tmpData.dataPtr(j)[i] << " ";
                osf << endl;
            }
            for (int i=0; i<nElts; ++i)
            {
                for (int j=0; j<MYLEN; ++j)
                    osf << connData[i*MYLEN+j] << " ";
                osf << endl;
            }
            osf << endl;
            osf.close();
        }
    }
    else
    {
        std::ofstream ofs;
        std::ostream* os =
            (outfile=="-" ? (std::ostream*)(&std::cout) : (std::ostream*)(&ofs) );
        if (outfile!="-")
            ofs.open(outfile.c_str(),std::ios::out|std::ios::trunc|std::ios::binary);
        (*os) << infile << " time = " << amrData.Time() << endl;
        (*os) << vars << endl;
        (*os) << nElts << " " << MYLEN << endl;
        tmpData.fab.writeOn(*os);
        (*os).write((char*)connData.dataPtr(),sizeof(int)*connData.size());
        if (outfile!="-")
            ofs.close();
    }

    BoxLib::Finalize();
    return 0;
}
Exemplo n.º 14
0
static
BoxArray
GetBndryCells (const BoxArray& ba,
               int             ngrow,
               const Geometry& geom)
{
    //
    // First get list of all ghost cells.
    //
    BoxList gcells, bcells;

    for (int i = 0; i < ba.size(); ++i)
	gcells.join(BoxLib::boxDiff(BoxLib::grow(ba[i],ngrow),ba[i]));
    //
    // Now strip out intersections with original BoxArray.
    //
    for (BoxList::const_iterator it = gcells.begin(); it != gcells.end(); ++it)
    {
        std::vector< std::pair<int,Box> > isects = ba.intersections(*it);

        if (isects.empty())
            bcells.push_back(*it);
        else
        {
            //
            // Collect all the intersection pieces.
            //
            BoxList pieces;
            for (int i = 0; i < isects.size(); i++)
                pieces.push_back(isects[i].second);
            BoxList leftover = BoxLib::complementIn(*it,pieces);
            bcells.catenate(leftover);
        }
    }
    //
    // Now strip out overlaps.
    //
    gcells.clear();
    gcells = BoxLib::removeOverlap(bcells);
    bcells.clear();

    if (geom.isAnyPeriodic())
    {
        Array<IntVect> pshifts(27);

        const Box& domain = geom.Domain();

        for (BoxList::const_iterator it = gcells.begin(); it != gcells.end(); ++it)
        {
            if (!domain.contains(*it))
            {
                //
                // Add in periodic ghost cells shifted to valid region.
                //
                geom.periodicShift(domain, *it, pshifts);

                for (int i = 0; i < pshifts.size(); i++)
                {
                    const Box& shftbox = *it + pshifts[i];

                    const Box& ovlp = domain & shftbox;
                    BoxList bl = BoxLib::complementIn(ovlp,BoxList(ba));
                    bcells.catenate(bl);
                }
            }
        }

        gcells.catenate(bcells);
    }

    return BoxArray(gcells);
}
Exemplo n.º 15
0
FabArrayBase::FBCacheIter
FabArrayBase::TheFB (bool                cross,
                     const FabArrayBase& mf)
{
    BL_PROFILE("FabArray::TheFB");

    BL_ASSERT(mf.size() > 0);

    const FabArrayBase::SI si(mf.boxArray(), mf.DistributionMap(), mf.nGrow(), cross);

    const IntVect& Typ   = mf.boxArray()[0].type();
    const int      Scale = D_TERM(Typ[0],+3*Typ[1],+5*Typ[2]) + 11;
    const int      Key   = mf.size() + mf.boxArray()[0].numPts() + mf.nGrow() + Scale + cross;

    std::pair<FBCacheIter,FBCacheIter> er_it = m_TheFBCache.equal_range(Key);

    for (FBCacheIter it = er_it.first; it != er_it.second; ++it)
    {
        if (it->second == si)
        {
	    ++it->second.m_nuse;
	    m_FBC_stats.recordUse();
            return it;
        }
    }

    if (m_TheFBCache.size() >= fb_cache_max_size)
    {
        //
        // Don't let the size of the cache get too big.
        // Get rid of entries with the biggest largest key that haven't been reused.
        // Otherwise just remove the entry with the largest key.
        //
        FBCacheIter End      = m_TheFBCache.end();
        FBCacheIter last_it  = End;
        FBCacheIter erase_it = End;

        for (FBCacheIter it = m_TheFBCache.begin(); it != End; ++it)
        {
            last_it = it;

            if (it->second.m_nuse <= 1)
                erase_it = it;
        }

        if (erase_it != End)
        {
	    m_FBC_stats.recordErase(erase_it->second.m_nuse);
            m_TheFBCache.erase(erase_it);
        }
        else if (last_it != End)
        {
	    m_FBC_stats.recordErase(last_it->second.m_nuse);
	    m_TheFBCache.erase(last_it);
        }
    }
    //
    // Got to insert one & then build it.
    //
    FBCacheIter                cache_it = m_TheFBCache.insert(FBCache::value_type(Key,si));
    SI&                        TheFB    = cache_it->second;
    const int                  MyProc   = ParallelDescriptor::MyProc();
    const BoxArray&            ba       = mf.boxArray();
    const DistributionMapping& dm       = mf.DistributionMap();
    const Array<int>&          imap     = mf.IndexMap();
    //
    // Here's where we allocate memory for the cache innards.
    // We do this so we don't have to build objects of these types
    // each time we search the cache.  Otherwise we'd be constructing
    // and destroying said objects quite frequently.
    //
    TheFB.m_LocTags = new CopyComTag::CopyComTagsContainer;
    TheFB.m_SndTags = new CopyComTag::MapOfCopyComTagContainers;
    TheFB.m_RcvTags = new CopyComTag::MapOfCopyComTagContainers;
    TheFB.m_SndVols = new std::map<int,int>;
    TheFB.m_RcvVols = new std::map<int,int>;

    TheFB.m_nuse = 1;

    m_FBC_stats.recordBuild();
    m_FBC_stats.recordUse();

    if (imap.empty())
        //
        // We don't own any of the relevant FABs so can't possibly have any work to do.
        //
        return cache_it;

    const int nlocal = imap.size();
    const int ng = si.m_ngrow;
    std::vector< std::pair<int,Box> > isects;

    CopyComTag::MapOfCopyComTagContainers send_tags; // temp copy

    for (int i = 0; i < nlocal; ++i)
    {
	const int ksnd = imap[i];
	const Box& vbx = ba[ksnd];

	ba.intersections(vbx, isects, ng);

	for (int j = 0, M = isects.size(); j < M; ++j)
	{
	    const int krcv      = isects[j].first;
	    const Box& bx       = isects[j].second;
	    const int dst_owner = dm[krcv];

	    if (krcv == ksnd) continue;  // same box

	    if (dst_owner == MyProc) continue;  // local copy will be dealt with later

	    send_tags[dst_owner].push_back(CopyComTag(bx, krcv, ksnd));
	}
    }

    CopyComTag::MapOfCopyComTagContainers recv_tags; // temp copy

    BaseFab<int> localtouch, remotetouch;
    bool check_local = false, check_remote = false;
#ifdef _OPENMP
    if (omp_get_max_threads() > 1) {
        check_local = true;
        check_remote = true;
    }
#endif

    if (ba.ixType().cellCentered()) {
	TheFB.m_threadsafe_loc = true;
	TheFB.m_threadsafe_rcv = true;
        check_local = false;
        check_remote = false;
    }

    for (int i = 0; i < nlocal; ++i)
    {
	const int   krcv = imap[i];
	const Box& bxrcv = BoxLib::grow(ba[krcv], ng);

	if (check_local) {
	    localtouch.resize(bxrcv);
	    localtouch.setVal(0);
	}

	if (check_remote) {
	    remotetouch.resize(bxrcv);
	    remotetouch.setVal(0);
	}

	ba.intersections(bxrcv, isects);

	for (int j = 0, M = isects.size(); j < M; ++j)
	{
	    const int ksnd      = isects[j].first;
	    const Box& bx       = isects[j].second;
	    const int src_owner = dm[ksnd];

	    if (krcv == ksnd) continue;  // same box

	    if (src_owner == MyProc) { // local copy
		const BoxList tilelist(bx, FabArrayBase::comm_tile_size);
		for (BoxList::const_iterator
			 it_tile  = tilelist.begin(),
			 End_tile = tilelist.end();   it_tile != End_tile; ++it_tile)
		{
		    TheFB.m_LocTags->push_back(CopyComTag(*it_tile, krcv, ksnd));
		}
		if (check_local) {
		    localtouch.plus(1, bx);
		}
	    } else {
		recv_tags[src_owner].push_back(CopyComTag(bx, krcv, ksnd));
		if (check_remote) {
		    remotetouch.plus(1, bx);
		}
	    }
	}

	if (check_local) {  
	    // safe if a cell is touched no more than once 
	    // keep checking thread safety if it is safe so far
            check_local = TheFB.m_threadsafe_loc = localtouch.max() <= 1;
        }

	if (check_remote) {
            check_remote = TheFB.m_threadsafe_rcv = remotetouch.max() <= 1;
        }
    }

//    ba.clear_hash_bin();

    for (int ipass = 0; ipass < 2; ++ipass) // pass 0: send; pass 1: recv
    {
	CopyComTag::MapOfCopyComTagContainers & Tags
	    = (ipass == 0) ? *TheFB.m_SndTags : *TheFB.m_RcvTags;
	CopyComTag::MapOfCopyComTagContainers & tmpTags
	    = (ipass == 0) ?        send_tags :        recv_tags;
	std::map<int,int> & Vols
	    = (ipass == 0) ? *TheFB.m_SndVols : *TheFB.m_RcvVols;

	for (CopyComTag::MapOfCopyComTagContainers::iterator 
		 it  = tmpTags.begin(), 
		 End = tmpTags.end();   it != End; ++it)
	{
	    const int key = it->first;
	    std::vector<CopyComTag>& cctv = it->second;

	    // We need to fix the order so that the send and recv processes match.
	    std::sort(cctv.begin(), cctv.end());

	    std::vector<CopyComTag> new_cctv;
	    new_cctv.reserve(cctv.size());

	    for (std::vector<CopyComTag>::const_iterator 
		     it2  = cctv.begin(),
		     End2 = cctv.end();   it2 != End2; ++it2)
	    {
		const Box& bx = it2->box;

		std::vector<Box> boxes;
		int vol = 0;

		if (si.m_cross) {
		    const Box& dstfabbx = ba[it2->fabIndex];
		    for (int dir = 0; dir < BL_SPACEDIM; dir++)
	            {
			Box lo = dstfabbx;
			lo.setSmall(dir, dstfabbx.smallEnd(dir) - ng);
			lo.setBig  (dir, dstfabbx.smallEnd(dir) - 1);
			lo &= bx;
			if (lo.ok()) {
			    boxes.push_back(lo);
			    vol += lo.numPts();
			}

			Box hi = dstfabbx;
			hi.setSmall(dir, dstfabbx.bigEnd(dir) + 1);
			hi.setBig  (dir, dstfabbx.bigEnd(dir) + ng);
			hi &= bx;
			if (hi.ok()) {
			    boxes.push_back(hi);
			    vol += hi.numPts();
			}
		    }
		} else {
		    boxes.push_back(bx);
		    vol += bx.numPts();
		}

		Vols[key] += vol;

		for (std::vector<Box>::const_iterator 
			 it_bx  = boxes.begin(),
			 End_bx = boxes.end();    it_bx != End_bx; ++it_bx)
	        {
		    const BoxList tilelist(*it_bx, FabArrayBase::comm_tile_size);
		    for (BoxList::const_iterator 
			     it_tile  = tilelist.begin(), 
			     End_tile = tilelist.end();   it_tile != End_tile; ++it_tile)
                    {
			new_cctv.push_back(CopyComTag(*it_tile, it2->fabIndex, it2->srcIndex));
		    }
		}
	    }

	    Tags[key].swap(new_cctv);
	}
    }

    return cache_it;
}
Exemplo n.º 16
0
void
MFGhostIter::Initialize ()
{
    int rit = 0;
    int nworkers = 1;
#ifdef BL_USE_TEAM
    if (ParallelDescriptor::TeamSize() > 1) {
	rit = ParallelDescriptor::MyRankInTeam();
	nworkers = ParallelDescriptor::TeamSize();
    }
#endif

    int tid = 0;
    int nthreads = 1;
#ifdef _OPENMP
    nthreads = omp_get_num_threads();
    if (nthreads > 1)
	tid = omp_get_thread_num();
#endif

    int npes = nworkers*nthreads;
    int pid = rit*nthreads+tid;

    BoxList alltiles;
    Array<int> allindex;
    Array<int> alllocalindex;

    for (int i=0; i < fabArray.IndexMap().size(); ++i) {
	int K = fabArray.IndexMap()[i];
	const Box& vbx = fabArray.box(K);
	const Box& fbx = fabArray.fabbox(K);

	const BoxList& diff = BoxLib::boxDiff(fbx, vbx);
	
	for (BoxList::const_iterator bli = diff.begin(); bli != diff.end(); ++bli) {
	    BoxList tiles(*bli, FabArrayBase::mfghostiter_tile_size);
	    int nt = tiles.size();
	    for (int it=0; it<nt; ++it) {
		allindex.push_back(K);
		alllocalindex.push_back(i);
	    }
	    alltiles.catenate(tiles);
	}
    }

    int n_tot_tiles = alltiles.size();
    int navg = n_tot_tiles / npes;
    int nleft = n_tot_tiles - navg*npes;
    int ntiles = navg;
    if (pid < nleft) ntiles++;

    // how many tiles should we skip?
    int nskip = pid*navg + std::min(pid,nleft);
    BoxList::const_iterator bli = alltiles.begin();
    for (int i=0; i<nskip; ++i) ++bli;

    lta.indexMap.reserve(ntiles);
    lta.localIndexMap.reserve(ntiles);
    lta.tileArray.reserve(ntiles);

    for (int i=0; i<ntiles; ++i) {
	lta.indexMap.push_back(allindex[i+nskip]);
	lta.localIndexMap.push_back(alllocalindex[i+nskip]);
	lta.tileArray.push_back(*bli++);
    }

    currentIndex = beginIndex = 0;
    endIndex = lta.indexMap.size();

    lta.nuse = 0;
    index_map       = &(lta.indexMap);
    local_index_map = &(lta.localIndexMap);
    tile_array      = &(lta.tileArray);
}