void MCLinOp::applyBC (MultiFab& inout, int level, MCBC_Mode bc_mode) { // // The inout MultiFab must have at least MCLinOp_grow ghost cells // for applyBC() // BL_ASSERT(inout.nGrow() >= MCLinOp_grow); // // The inout MultiFab must have at least Periodic_BC_grow cells for the // algorithms taking care of periodic boundary conditions. // BL_ASSERT(inout.nGrow() >= MCLinOp_grow); // // No coarsened boundary values, cannot apply inhomog at lev>0. // BL_ASSERT(!(level>0 && bc_mode == MCInhomogeneous_BC)); int flagden = 1; // fill in the bndry data and undrrelxr int flagbc = 1; // with values if (bc_mode == MCHomogeneous_BC) flagbc = 0; // nodata if homog int nc = inout.nComp(); BL_ASSERT(nc == numcomp ); inout.setBndry(-1.e30); inout.FillBoundary(); prepareForLevel(level); geomarray[level].FillPeriodicBoundary(inout,0,nc); // // Fill boundary cells. // #ifdef _OPENMP #pragma omp parallel #endif for (MFIter mfi(inout); mfi.isValid(); ++mfi) { const int gn = mfi.index(); BL_ASSERT(gbox[level][gn] == inout.box(gn)); const BndryData::RealTuple& bdl = bgb.bndryLocs(gn); const Array< Array<BoundCond> >& bdc = bgb.bndryConds(gn); const MaskTuple& msk = maskvals[level][gn]; for (OrientationIter oitr; oitr; ++oitr) { const Orientation face = oitr(); FabSet& f = (*undrrelxr[level])[face]; FabSet& td = (*tangderiv[level])[face]; int cdr(face); const FabSet& fs = bgb.bndryValues(face); Real bcl = bdl[face]; const Array<BoundCond>& bc = bdc[face]; const int *bct = (const int*) bc.dataPtr(); const FArrayBox& fsfab = fs[gn]; const Real* bcvalptr = fsfab.dataPtr(); // // Way external derivs stored. // const Real* exttdptr = fsfab.dataPtr(numcomp); const int* fslo = fsfab.loVect(); const int* fshi = fsfab.hiVect(); FArrayBox& inoutfab = inout[gn]; FArrayBox& denfab = f[gn]; FArrayBox& tdfab = td[gn]; #if BL_SPACEDIM==2 int cdir = face.coordDir(), perpdir = -1; if (cdir == 0) perpdir = 1; else if (cdir == 1) perpdir = 0; else BoxLib::Abort("MCLinOp::applyBC(): bad logic"); const Mask& m = *msk[face]; const Mask& mphi = *msk[Orientation(perpdir,Orientation::high)]; const Mask& mplo = *msk[Orientation(perpdir,Orientation::low)]; FORT_APPLYBC( &flagden, &flagbc, &maxorder, inoutfab.dataPtr(), ARLIM(inoutfab.loVect()), ARLIM(inoutfab.hiVect()), &cdr, bct, &bcl, bcvalptr, ARLIM(fslo), ARLIM(fshi), m.dataPtr(), ARLIM(m.loVect()), ARLIM(m.hiVect()), mphi.dataPtr(), ARLIM(mphi.loVect()), ARLIM(mphi.hiVect()), mplo.dataPtr(), ARLIM(mplo.loVect()), ARLIM(mplo.hiVect()), denfab.dataPtr(), ARLIM(denfab.loVect()), ARLIM(denfab.hiVect()), exttdptr, ARLIM(fslo), ARLIM(fshi), tdfab.dataPtr(),ARLIM(tdfab.loVect()),ARLIM(tdfab.hiVect()), inout.box(gn).loVect(), inout.box(gn).hiVect(), &nc, h[level]); #elif BL_SPACEDIM==3 const Mask& mn = *msk[Orientation(1,Orientation::high)]; const Mask& me = *msk[Orientation(0,Orientation::high)]; const Mask& mw = *msk[Orientation(0,Orientation::low)]; const Mask& ms = *msk[Orientation(1,Orientation::low)]; const Mask& mt = *msk[Orientation(2,Orientation::high)]; const Mask& mb = *msk[Orientation(2,Orientation::low)]; FORT_APPLYBC( &flagden, &flagbc, &maxorder, inoutfab.dataPtr(), ARLIM(inoutfab.loVect()), ARLIM(inoutfab.hiVect()), &cdr, bct, &bcl, bcvalptr, ARLIM(fslo), ARLIM(fshi), mn.dataPtr(),ARLIM(mn.loVect()),ARLIM(mn.hiVect()), me.dataPtr(),ARLIM(me.loVect()),ARLIM(me.hiVect()), mw.dataPtr(),ARLIM(mw.loVect()),ARLIM(mw.hiVect()), ms.dataPtr(),ARLIM(ms.loVect()),ARLIM(ms.hiVect()), mt.dataPtr(),ARLIM(mt.loVect()),ARLIM(mt.hiVect()), mb.dataPtr(),ARLIM(mb.loVect()),ARLIM(mb.hiVect()), denfab.dataPtr(), ARLIM(denfab.loVect()), ARLIM(denfab.hiVect()), exttdptr, ARLIM(fslo), ARLIM(fshi), tdfab.dataPtr(),ARLIM(tdfab.loVect()),ARLIM(tdfab.hiVect()), inout.box(gn).loVect(), inout.box(gn).hiVect(), &nc, h[level]); #endif } } #if 0 // This "probably" works, but is not strictly needed just because of the way Bill // coded up the tangential derivative stuff. It's handy code though, so I want to // keep it around/ // Clean up corners: // The problem here is that APPLYBC fills only grow cells normal to the boundary. // As a result, any corner cell on the boundary (either coarse-fine or fine-fine) // is not filled. For coarse-fine, the operator adjusts itself, sliding away from // the box edge to avoid referencing that corner point. On the physical boundary // though, the corner point is needed. Particularly if a fine-fine boundary intersects // the physical boundary, since we want the stencil to be independent of the box // blocking. FillBoundary operations wont fix the problem because the "good" // data we need is living in the grow region of adjacent fabs. So, here we play // the usual games to treat the newly filled grow cells as "valid" data. // Note that we only need to do something where the grids touch the physical boundary. const Geometry& geomlev = geomarray[level]; const BoxArray& grids = inout.boxArray(); const Box& domain = geomlev.Domain(); int nGrow = 1; int src_comp = 0; int num_comp = BL_SPACEDIM; // Lets do a quick check to see if we need to do anything at all here BoxArray BIGba = BoxArray(grids).grow(nGrow); if (! (domain.contains(BIGba.minimalBox())) ) { BoxArray boundary_pieces; Array<int> proc_idxs; Array<Array<int> > old_to_new(grids.size()); const DistributionMapping& dmap=inout.DistributionMap(); for (int d=0; d<BL_SPACEDIM; ++d) { if (! (geomlev.isPeriodic(d)) ) { BoxArray gba = BoxArray(grids).grow(d,nGrow); for (int i=0; i<gba.size(); ++i) { BoxArray new_pieces = BoxLib::boxComplement(gba[i],domain); int size_new = new_pieces.size(); if (size_new>0) { int size_old = boundary_pieces.size(); boundary_pieces.resize(size_old+size_new); proc_idxs.resize(boundary_pieces.size()); for (int j=0; j<size_new; ++j) { boundary_pieces.set(size_old+j,new_pieces[j]); proc_idxs[size_old+j] = dmap[i]; old_to_new[i].push_back(size_old+j); } } } } } proc_idxs.push_back(ParallelDescriptor::MyProc()); MultiFab boundary_data(boundary_pieces,num_comp,nGrow, DistributionMapping(proc_idxs)); for (MFIter mfi(inout); mfi.isValid(); ++mfi) { const FArrayBox& src_fab = inout[mfi]; for (int j=0; j<old_to_new[mfi.index()].size(); ++j) { int new_box_idx = old_to_new[mfi.index()][j]; boundary_data[new_box_idx].copy(src_fab,src_comp,0,num_comp); } } boundary_data.FillBoundary(); // Use a hacked Geometry object to handle the periodic intersections for us. // Here, the "domain" is the plane of cells on non-periodic boundary faces. // and there may be cells over the periodic boundary in the remaining directions. // We do a Geometry::PFB on each non-periodic face to sync these up. if (geomlev.isAnyPeriodic()) { Array<int> is_per(BL_SPACEDIM,0); for (int d=0; d<BL_SPACEDIM; ++d) { is_per[d] = geomlev.isPeriodic(d); } for (int d=0; d<BL_SPACEDIM; ++d) { if (! is_per[d]) { Box tmpLo = BoxLib::adjCellLo(geomlev.Domain(),d,1); Geometry tmpGeomLo(tmpLo,&(geomlev.ProbDomain()),(int)geomlev.Coord(),is_per.dataPtr()); tmpGeomLo.FillPeriodicBoundary(boundary_data); Box tmpHi = BoxLib::adjCellHi(geomlev.Domain(),d,1); Geometry tmpGeomHi(tmpHi,&(geomlev.ProbDomain()),(int)geomlev.Coord(),is_per.dataPtr()); tmpGeomHi.FillPeriodicBoundary(boundary_data); } } } for (MFIter mfi(inout); mfi.isValid(); ++mfi) { int idx = mfi.index(); FArrayBox& dst_fab = inout[mfi]; for (int j=0; j<old_to_new[idx].size(); ++j) { int new_box_idx = old_to_new[mfi.index()][j]; const FArrayBox& src_fab = boundary_data[new_box_idx]; const Box& src_box = src_fab.box(); BoxArray pieces_outside_domain = BoxLib::boxComplement(src_box,domain); for (int k=0; k<pieces_outside_domain.size(); ++k) { const Box& outside = pieces_outside_domain[k] & dst_fab.box(); if (outside.ok()) { dst_fab.copy(src_fab,outside,0,outside,src_comp,num_comp); } } } } } #endif }
int main(int argc, char* argv[]) { BoxLib::Initialize(argc,argv); BL_PROFILE_VAR("main()", pmain); std::cout << std::setprecision(15); ParmParse ppmg("mg"); ppmg.query("v", verbose); ppmg.query("maxorder", maxorder); ParmParse pp; { std::string solver_type_s; pp.get("solver_type",solver_type_s); if (solver_type_s == "BoxLib_C") { solver_type = BoxLib_C; } else if (solver_type_s == "BoxLib_C4") { solver_type = BoxLib_C4; } else if (solver_type_s == "BoxLib_F") { #ifdef USE_F90_SOLVERS solver_type = BoxLib_F; #else BoxLib::Error("Set USE_FORTRAN=TRUE in GNUmakefile"); #endif } else if (solver_type_s == "Hypre") { #ifdef USEHYPRE solver_type = Hypre; #else BoxLib::Error("Set USE_HYPRE=TRUE in GNUmakefile"); #endif } else if (solver_type_s == "All") { solver_type = All; } else { if (ParallelDescriptor::IOProcessor()) { std::cout << "Don't know this solver type: " << solver_type << std::endl; } BoxLib::Error(""); } } { std::string bc_type_s; pp.get("bc_type",bc_type_s); if (bc_type_s == "Dirichlet") { bc_type = Dirichlet; #ifdef USEHPGMG domain_boundary_condition = BC_DIRICHLET; #endif } else if (bc_type_s == "Neumann") { bc_type = Neumann; #ifdef USEHPGMG BoxLib::Error("HPGMG does not support Neumann boundary conditions"); #endif } else if (bc_type_s == "Periodic") { bc_type = Periodic; #ifdef USEHPGMG domain_boundary_condition = BC_PERIODIC; #endif } else { if (ParallelDescriptor::IOProcessor()) { std::cout << "Don't know this boundary type: " << bc_type << std::endl; } BoxLib::Error(""); } } pp.query("tol_rel", tolerance_rel); pp.query("tol_abs", tolerance_abs); pp.query("maxiter", maxiter); pp.query("plot_rhs" , plot_rhs); pp.query("plot_beta", plot_beta); pp.query("plot_soln", plot_soln); pp.query("plot_asol", plot_asol); pp.query("plot_err", plot_err); pp.query("comp_norm", comp_norm); Real a, b; pp.get("a", a); pp.get("b", b); pp.get("n_cell",n_cell); pp.get("max_grid_size",max_grid_size); // Define a single box covering the domain IntVect dom_lo(D_DECL(0,0,0)); IntVect dom_hi(D_DECL(n_cell-1,n_cell-1,n_cell-1)); Box domain(dom_lo,dom_hi); // Initialize the boxarray "bs" from the single box "bx" BoxArray bs(domain); // Break up boxarray "bs" into chunks no larger than "max_grid_size" along a direction bs.maxSize(max_grid_size); // This defines the physical size of the box. Right now the box is [0,1] in each direction. RealBox real_box; for (int n = 0; n < BL_SPACEDIM; n++) { real_box.setLo(n, 0.0); real_box.setHi(n, 1.0); } // This says we are using Cartesian coordinates int coord = 0; // This sets the boundary conditions to be periodic or not Array<int> is_per(BL_SPACEDIM,1); if (bc_type == Dirichlet || bc_type == Neumann) { if (ParallelDescriptor::IOProcessor()) { std::cout << "Using Dirichlet or Neumann boundary conditions." << std::endl; } for (int n = 0; n < BL_SPACEDIM; n++) is_per[n] = 0; } else { if (ParallelDescriptor::IOProcessor()) { std::cout << "Using periodic boundary conditions." << std::endl; } for (int n = 0; n < BL_SPACEDIM; n++) is_per[n] = 1; } // This defines a Geometry object which is useful for writing the plotfiles Geometry geom(domain,&real_box,coord,is_per.dataPtr()); for ( int n=0; n<BL_SPACEDIM; n++ ) { dx[n] = ( geom.ProbHi(n) - geom.ProbLo(n) )/domain.length(n); } if (ParallelDescriptor::IOProcessor()) { std::cout << "Grid resolution : " << n_cell << " (cells)" << std::endl; std::cout << "Domain size : " << real_box.hi(0) - real_box.lo(0) << " (length unit) " << std::endl; std::cout << "Max_grid_size : " << max_grid_size << " (cells)" << std::endl; std::cout << "Number of grids : " << bs.size() << std::endl; } // Allocate and define the right hand side. bool do_4th = (solver_type==BoxLib_C4 || solver_type==All); int ngr = (do_4th ? 1 : 0); MultiFab rhs(bs, Ncomp, ngr); setup_rhs(rhs, geom); // Set up the Helmholtz operator coefficients. MultiFab alpha(bs, Ncomp, 0); PArray<MultiFab> beta(BL_SPACEDIM, PArrayManage); for ( int n=0; n<BL_SPACEDIM; ++n ) { BoxArray bx(bs); beta.set(n, new MultiFab(bx.surroundingNodes(n), Ncomp, 0, Fab_allocate)); } // The way HPGMG stores face-centered data is completely different than the // way BoxLib does it, and translating between the two directly via indexing // magic is a nightmare. Happily, the way this tutorial calculates // face-centered values is by first calculating cell-centered values and then // interpolating to the cell faces. HPGMG can do the same thing, so rather // than converting directly from BoxLib's face-centered data to HPGMG's, just // give HPGMG the cell-centered data and let it interpolate itself. MultiFab beta_cc(bs,Ncomp,1); // cell-centered beta setup_coeffs(bs, alpha, beta, geom, beta_cc); MultiFab alpha4, beta4; if (do_4th) { alpha4.define(bs, Ncomp, 4, Fab_allocate); beta4.define(bs, Ncomp, 3, Fab_allocate); setup_coeffs4(bs, alpha4, beta4, geom); } MultiFab anaSoln; if (comp_norm || plot_err || plot_asol) { anaSoln.define(bs, Ncomp, 0, Fab_allocate); compute_analyticSolution(anaSoln,Array<Real>(BL_SPACEDIM,0.5)); if (plot_asol) { writePlotFile("ASOL", anaSoln, geom); } } // Allocate the solution array // Set the number of ghost cells in the solution array. MultiFab soln(bs, Ncomp, 1); MultiFab soln4; if (do_4th) { soln4.define(bs, Ncomp, 3, Fab_allocate); } MultiFab gphi(bs, BL_SPACEDIM, 0); #ifdef USEHYPRE if (solver_type == Hypre || solver_type == All) { if (ParallelDescriptor::IOProcessor()) { std::cout << "----------------------------------------" << std::endl; std::cout << "Solving with Hypre " << std::endl; } solve(soln, anaSoln, gphi, a, b, alpha, beta, beta_cc, rhs, bs, geom, Hypre); } #endif if (solver_type == BoxLib_C || solver_type == All) { if (ParallelDescriptor::IOProcessor()) { std::cout << "----------------------------------------" << std::endl; std::cout << "Solving with BoxLib C++ solver " << std::endl; } solve(soln, anaSoln, gphi, a, b, alpha, beta, beta_cc, rhs, bs, geom, BoxLib_C); } if (solver_type == BoxLib_C4 || solver_type == All) { if (ParallelDescriptor::IOProcessor()) { std::cout << "----------------------------------------" << std::endl; std::cout << "Solving with BoxLib C++ 4th order solver " << std::endl; } solve4(soln4, anaSoln, a, b, alpha4, beta4, rhs, bs, geom); } #ifdef USE_F90_SOLVERS if (solver_type == BoxLib_F || solver_type == All) { if (ParallelDescriptor::IOProcessor()) { std::cout << "----------------------------------------" << std::endl; std::cout << "Solving with BoxLib F90 solver " << std::endl; } solve(soln, anaSoln, gphi, a, b, alpha, beta, beta_cc, rhs, bs, geom, BoxLib_F); } #endif #ifdef USEHPGMG if (solver_type == HPGMG || solver_type == All) { if (ParallelDescriptor::IOProcessor()) { std::cout << "----------------------------------------" << std::endl; std::cout << "Solving with HPGMG solver " << std::endl; } solve(soln, anaSoln, gphi, a, b, alpha, beta, beta_cc, rhs, bs, geom, HPGMG); } #endif if (ParallelDescriptor::IOProcessor()) { std::cout << "----------------------------------------" << std::endl; } BL_PROFILE_VAR_STOP(pmain); BoxLib::Finalize(); }