Esempio n. 1
0
void
MCLinOp::makeCoefficients (MultiFab&       cs,
                           const MultiFab& fn,
                           int             level)
{
    const int nc = fn.nComp();
    //
    // Determine index type of incoming MultiFab.
    //
    const IndexType iType(fn.boxArray().ixType());
    const IndexType cType(D_DECL(IndexType::CELL, IndexType::CELL, IndexType::CELL));
    const IndexType xType(D_DECL(IndexType::NODE, IndexType::CELL, IndexType::CELL));
    const IndexType yType(D_DECL(IndexType::CELL, IndexType::NODE, IndexType::CELL));
#if (BL_SPACEDIM == 3)    
    const IndexType zType(D_DECL(IndexType::CELL, IndexType::CELL, IndexType::NODE));
#endif
    int cdir;
    if (iType == cType)
    {
        cdir = -1;
    }
    else if (iType == xType)
    {
        cdir = 0;
    }
    else if (iType == yType)
    {
        cdir = 1;
    }
#if (BL_SPACEDIM == 3)
    else if (iType == zType)
    {
        cdir = 2;
    }
#endif
    else
        BoxLib::Abort("MCLinOp::makeCoeffients(): Bad index type");
    
    BoxArray d(gbox[level]);
    if (cdir >= 0)
	d.surroundingNodes(cdir);

    int nGrow=0;
    cs.define(d, nc, nGrow, Fab_allocate);
    cs.setVal(0.0);

    const BoxArray& grids = gbox[level];

    for (MFIter csmfi(cs); csmfi.isValid(); ++csmfi)
    {
        const Box&       grd   = grids[csmfi.index()];
        FArrayBox&       csfab = cs[csmfi];
        const FArrayBox& fnfab = fn[csmfi];

	switch(cdir)
        {
	case -1:
	    FORT_AVERAGECC(
		csfab.dataPtr(),
                ARLIM(csfab.loVect()), ARLIM(csfab.hiVect()),
		fnfab.dataPtr(),
                ARLIM(fnfab.loVect()), ARLIM(fnfab.hiVect()),
		grd.loVect(),
                grd.hiVect(), &nc);
	    break;
	case 0:
	case 1:
	case 2:
	    if ( harmavg )
            {
		FORT_HARMONIC_AVERAGEEC(
		    csfab.dataPtr(), 
                    ARLIM(csfab.loVect()), ARLIM(csfab.hiVect()),
		    fnfab.dataPtr(), 
                    ARLIM(fnfab.loVect()), ARLIM(fnfab.hiVect()),
		    grd.loVect(),
                    grd.hiVect(), &nc, &cdir);
	    }
            else
            {
		FORT_AVERAGEEC(
		    csfab.dataPtr(), 
                    ARLIM(csfab.loVect()), ARLIM(csfab.hiVect()),
		    fnfab.dataPtr(), 
                    ARLIM(fnfab.loVect()), ARLIM(fnfab.hiVect()),
		    grd.loVect(),
                    grd.hiVect(), &nc, &cdir);
	    }
	    break;
	default:
	    BoxLib::Error("MCLinOp::makeCoeffients(): bad coefficient coarsening direction!");
	}
    }
}
Esempio n. 2
0
int main(int argc, char* argv[])
{
  BoxLib::Initialize(argc,argv);

  BL_PROFILE_VAR("main()", pmain);

  std::cout << std::setprecision(15);

  solver_type = BoxLib_C;
  bc_type = Periodic;

  Real     a = 0.0;
  Real     b = 1.0;

  // ---- First use the number of processors to decide how many grids you have.
  // ---- We arbitrarily decide to have one grid per MPI process in a uniform
  // ---- cubic domain, so we require that the number of processors be N^3.
  // ---- This requirement is somewhat arbitrary, but convenient for now.

  int nprocs = ParallelDescriptor::NProcs();

  // N is the cube root of the number of processors
  int N(0);
  for(int i(1); i*i*i <= nprocs; ++i) {
    if(i*i*i == nprocs) {
      N = i;
    }
  }

  if(N == 0) {  // not a cube
    if(ParallelDescriptor::IOProcessor()) {
      std::cerr << "**** Error:  nprocs = " << nprocs << " is not currently supported." << std::endl;
    }
    BoxLib::Error("We require that the number of processors be a perfect cube");
  }
  if(ParallelDescriptor::IOProcessor()) {
    std::cout << "N = " << N << std::endl;
  }


  // ---- make a box, then a boxarray with maxSize
  int domain_hi = (N*maxGrid) - 1;
  Box domain(IntVect(0,0,0), IntVect(domain_hi,domain_hi,domain_hi));
  BoxArray bs(domain);
  bs.maxSize(maxGrid);

  // 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
  int is_per[BL_SPACEDIM];
  
  if (bc_type == Dirichlet || bc_type == Neumann) {
    for (int n = 0; n < BL_SPACEDIM; n++) is_per[n] = 0;
  } 
  else {
    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);

  for ( int n=0; n<BL_SPACEDIM; n++ ) {
    dx[n] = ( geom.ProbHi(n) - geom.ProbLo(n) )/domain.length(n);
  }

  if (ParallelDescriptor::IOProcessor()) {
     std::cout << "Domain size     : " << N << std::endl;
     std::cout << "Max_grid_size   : " << maxGrid << std::endl;
     std::cout << "Number of grids : " << bs.size() << std::endl;
  }

  // Allocate and define the right hand side.
  MultiFab rhs(bs, Ncomp, 0, Fab_allocate); 
  setup_rhs(rhs, geom, a, b);

  MultiFab alpha(bs, Ncomp, 0, Fab_allocate);
  MultiFab beta[BL_SPACEDIM];
  for ( int n=0; n<BL_SPACEDIM; ++n ) {
    BoxArray bx(bs);
    beta[n].define(bx.surroundingNodes(n), Ncomp, 1, Fab_allocate);
  }

  setup_coeffs(bs, alpha, beta, geom);

  MultiFab anaSoln;
  if (comp_norm) {
    anaSoln.define(bs, Ncomp, 0, Fab_allocate);
    compute_analyticSolution(anaSoln);
  }

  // Allocate the solution array 
  // Set the number of ghost cells in the solution array.
  MultiFab soln(bs, Ncomp, 1, Fab_allocate);

  solve(soln, anaSoln, a, b, alpha, beta, rhs, bs, geom, BoxLib_C);

  BL_PROFILE_VAR_STOP(pmain);

  BoxLib::Finalize();
}
Esempio n. 3
0
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();
}
Esempio n. 4
0
void
LinOp::makeCoefficients (MultiFab&       cs,
                         const MultiFab& fn,
                         int             level)
{
    BL_PROFILE("LinOp::makeCoefficients()");

    int nc = 1;
    //
    // Determine index type of incoming MultiFab.
    //
    const IndexType iType(fn.boxArray().ixType());
    const IndexType cType(D_DECL(IndexType::CELL, IndexType::CELL, IndexType::CELL));
    const IndexType xType(D_DECL(IndexType::NODE, IndexType::CELL, IndexType::CELL));
    const IndexType yType(D_DECL(IndexType::CELL, IndexType::NODE, IndexType::CELL));
#if (BL_SPACEDIM == 3)    
    const IndexType zType(D_DECL(IndexType::CELL, IndexType::CELL, IndexType::NODE));
#endif

    int cdir;
    if (iType == cType)
    {
        cdir = -1;
    }
    else if (iType == xType)
    {
        cdir = 0;
    }
    else if (iType == yType)
    {
        cdir = 1;
#if (BL_SPACEDIM == 3)
    }
    else if (iType == zType)
    {
        cdir = 2;
#endif    
    }
    else
    {
        BoxLib::Error("LinOp::makeCoeffients: Bad index type");
    }

    BoxArray d(gbox[level]);
    if (cdir >= 0)
        d.surroundingNodes(cdir);
    //
    // Only single-component solves supported (verified) by this class.
    //
    const int nComp=1;
    const int nGrow=0;
    cs.define(d, nComp, nGrow, Fab_allocate);

    const bool tiling = true;

    switch (cdir)
    {
    case -1:
#ifdef _OPENMP
#pragma omp parallel
#endif
        for (MFIter csmfi(cs,tiling); csmfi.isValid(); ++csmfi)
        {
            const Box& tbx = csmfi.tilebox();
            FArrayBox&       csfab = cs[csmfi];
            const FArrayBox& fnfab = fn[csmfi];

            FORT_AVERAGECC(csfab.dataPtr(), ARLIM(csfab.loVect()),
                           ARLIM(csfab.hiVect()),fnfab.dataPtr(),
                           ARLIM(fnfab.loVect()),ARLIM(fnfab.hiVect()),
                           tbx.loVect(),tbx.hiVect(), &nc);
        }
        break;
    case 0:
    case 1:
    case 2:
        if (harmavg)
        {
#ifdef _OPENMP
#pragma omp parallel
#endif
  	    for (MFIter csmfi(cs,tiling); csmfi.isValid(); ++csmfi)
            {
	        const Box& tbx = csmfi.tilebox();
                FArrayBox&       csfab = cs[csmfi];
                const FArrayBox& fnfab = fn[csmfi];

                FORT_HARMONIC_AVERAGEEC(csfab.dataPtr(),
                                        ARLIM(csfab.loVect()),
                                        ARLIM(csfab.hiVect()),
                                        fnfab.dataPtr(),
                                        ARLIM(fnfab.loVect()),
                                        ARLIM(fnfab.hiVect()),
                                        tbx.loVect(),tbx.hiVect(),
                                        &nc,&cdir);
            }
        }
        else
        {
#ifdef _OPENMP
#pragma omp parallel
#endif
            for (MFIter csmfi(cs,tiling); csmfi.isValid(); ++csmfi)
            {
                const Box& tbx = csmfi.tilebox();
                FArrayBox&       csfab = cs[csmfi];
                const FArrayBox& fnfab = fn[csmfi];

                FORT_AVERAGEEC(csfab.dataPtr(),ARLIM(csfab.loVect()),
                               ARLIM(csfab.hiVect()),fnfab.dataPtr(), 
                               ARLIM(fnfab.loVect()),ARLIM(fnfab.hiVect()),
	                       tbx.loVect(),tbx.hiVect(),
                               &nc, &cdir);
            }
        }
        break;
    default:
        BoxLib::Error("LinOp:: bad coefficient coarsening direction!");
    }
}
Esempio n. 5
0
int
main (int   argc,
      char* argv[])
{
    BoxLib::Initialize(argc,argv);

    std::cout << std::setprecision(10);

    if (argc < 2)
    {
      std::cerr << "usage:  " << argv[0] << " inputsfile [options]" << '\n';
      exit(-1);
    }

    ParmParse pp;
    
    int n;

    BoxArray bs;
    
#if BL_SPACEDIM == 2
    Box domain(IntVect(0,0),IntVect(11,11));
    std::string boxfile("gr.2_small_a") ;
#elif BL_SPACEDIM == 3
    Box domain(IntVect(0,0,0),IntVect(11,11,11));
    std::string boxfile("grids/gr.3_2x3x4") ;
#endif
    pp.query("boxes", boxfile);

    std::ifstream ifs(boxfile.c_str(), std::ios::in);

    if (!ifs)
    {
        std::string msg = "problem opening grids file: ";
        msg += boxfile.c_str();
        BoxLib::Abort(msg.c_str());
    }

    ifs >> domain;

    if (ParallelDescriptor::IOProcessor())
	std::cout << "domain: " << domain << std::endl;

    bs.readFrom(ifs);

    if (ParallelDescriptor::IOProcessor())
	std::cout << "grids:\n" << bs << std::endl;

    Geometry geom(domain);
    const Real* H = geom.CellSize();
    int ratio=2; pp.query("ratio", ratio);

    // allocate/init soln and rhs
    int Ncomp=BL_SPACEDIM;
    int Nghost=0;
    int Ngrids=bs.size();
    MultiFab soln(bs, Ncomp, Nghost, Fab_allocate); soln.setVal(0.0);
    MultiFab out(bs, Ncomp, Nghost, Fab_allocate); 
    MultiFab rhs(bs, Ncomp, Nghost, Fab_allocate); rhs.setVal(0.0);
    for(MFIter rhsmfi(rhs); rhsmfi.isValid(); ++rhsmfi)
    {
	FORT_FILLRHS(rhs[rhsmfi].dataPtr(),
		     ARLIM(rhs[rhsmfi].loVect()),ARLIM(rhs[rhsmfi].hiVect()),
		     H,&Ncomp);
    }
    
    // Create the boundary object
    MCViscBndry vbd(bs,geom);

    BCRec phys_bc;
    Array<int> lo_bc(BL_SPACEDIM), hi_bc(BL_SPACEDIM);
    pp.getarr("lo_bc",lo_bc,0,BL_SPACEDIM);
    pp.getarr("hi_bc",hi_bc,0,BL_SPACEDIM);
    for (int i = 0; i < BL_SPACEDIM; i++)
    {
        phys_bc.setLo(i,lo_bc[i]);
        phys_bc.setHi(i,hi_bc[i]);
    }

    
    // Create the BCRec's interpreted by ViscBndry objects
#if BL_SPACEDIM==2
    Array<BCRec> pbcarray(4);
    pbcarray[0] = BCRec(D_DECL(REFLECT_ODD,REFLECT_EVEN,EXT_DIR),
			D_DECL(EXT_DIR,EXT_DIR,EXT_DIR));
    pbcarray[1] = BCRec(D_DECL(REFLECT_EVEN,REFLECT_ODD,EXT_DIR),
			D_DECL(EXT_DIR,EXT_DIR,EXT_DIR));
    pbcarray[2] = BCRec(D_DECL(EXT_DIR,EXT_DIR,EXT_DIR),
			D_DECL(EXT_DIR,EXT_DIR,EXT_DIR));
    pbcarray[3] = BCRec(D_DECL(EXT_DIR,EXT_DIR,EXT_DIR),
			D_DECL(EXT_DIR,EXT_DIR,EXT_DIR));
#elif BL_SPACEDIM==3
    Array<BCRec> pbcarray(12);

#if 1
    pbcarray[0] = BCRec(EXT_DIR,EXT_DIR,EXT_DIR,EXT_DIR,EXT_DIR,EXT_DIR);
    pbcarray[1] = BCRec(EXT_DIR,EXT_DIR,EXT_DIR,EXT_DIR,EXT_DIR,EXT_DIR);
    pbcarray[2] = BCRec(EXT_DIR,EXT_DIR,EXT_DIR,EXT_DIR,EXT_DIR,EXT_DIR);
    pbcarray[3] = BCRec(D_DECL(EXT_DIR,EXT_DIR,EXT_DIR),
			D_DECL(EXT_DIR,EXT_DIR,EXT_DIR));
    pbcarray[4] = BCRec(D_DECL(EXT_DIR,EXT_DIR,EXT_DIR),
			D_DECL(EXT_DIR,EXT_DIR,EXT_DIR));
    pbcarray[5] = BCRec(D_DECL(EXT_DIR,EXT_DIR,EXT_DIR),
			D_DECL(EXT_DIR,EXT_DIR,EXT_DIR));
    pbcarray[6] = BCRec(D_DECL(EXT_DIR,EXT_DIR,EXT_DIR),
			D_DECL(EXT_DIR,EXT_DIR,EXT_DIR));
    pbcarray[7] = BCRec(D_DECL(EXT_DIR,EXT_DIR,EXT_DIR),
			D_DECL(EXT_DIR,EXT_DIR,EXT_DIR));
    pbcarray[8] = BCRec(D_DECL(EXT_DIR,EXT_DIR,EXT_DIR),
			D_DECL(EXT_DIR,EXT_DIR,EXT_DIR));
    pbcarray[9] = BCRec(D_DECL(EXT_DIR,EXT_DIR,EXT_DIR),
			D_DECL(EXT_DIR,EXT_DIR,EXT_DIR));
    pbcarray[10] = BCRec(D_DECL(EXT_DIR,EXT_DIR,EXT_DIR),
			 D_DECL(EXT_DIR,EXT_DIR,EXT_DIR));
    pbcarray[11] = BCRec(D_DECL(EXT_DIR,EXT_DIR,EXT_DIR),
			 D_DECL(EXT_DIR,EXT_DIR,EXT_DIR));
#else
    for (int i = 0; i < 12; i++)
        pbcarray[i] = phys_bc;
#endif
#endif
    
    Nghost = 1; // need space for bc info
    MultiFab fine(bs,Ncomp,Nghost,Fab_allocate);
    for(MFIter finemfi(fine); finemfi.isValid(); ++finemfi)
    {
	FORT_FILLFINE(fine[finemfi].dataPtr(),
		      ARLIM(fine[finemfi].loVect()),ARLIM(fine[finemfi].hiVect()),
		      H,&Ncomp);
    }

    // Create "background coarse data"
    Box crse_bx = Box(domain).coarsen(ratio).grow(1);
    BoxArray cba(crse_bx);
    cba.maxSize(32);
    Real h_crse[BL_SPACEDIM];
    for (n=0; n<BL_SPACEDIM; n++) h_crse[n] = H[n]*ratio;

    MultiFab crse_mf(cba, Ncomp, 0);
//    FArrayBox crse_fab(crse_bx,Ncomp);

    for (MFIter mfi(crse_mf); mfi.isValid(); ++mfi)
    {
        FORT_FILLCRSE(crse_mf[mfi].dataPtr(),
                      ARLIM(crse_mf[mfi].loVect()),ARLIM(crse_mf[mfi].hiVect()),
                      h_crse,&Ncomp);
    }


    
    // Create coarse boundary register, fill w/data from coarse FAB
    int bndry_InRad=0;
    int bndry_OutRad=1;
    int bndry_Extent=1;
    BoxArray cbs = BoxArray(bs).coarsen(ratio);
    BndryRegister cbr(cbs,bndry_InRad,bndry_OutRad,bndry_Extent,Ncomp);
    for (OrientationIter face; face; ++face)
    {
	Orientation f = face();
	FabSet& bnd_fs(cbr[f]);
	bnd_fs.copyFrom(crse_mf, 0, 0, 0, Ncomp);
    }
  
    // Interpolate crse data to fine boundary, where applicable
    int cbr_Nstart=0;
    int fine_Nstart=0;
    int bndry_Nstart=0;
    vbd.setBndryValues(cbr,cbr_Nstart,fine,fine_Nstart,
		       bndry_Nstart,Ncomp,ratio,pbcarray);
  
    Nghost = 1; // other variables don't need extra space
    
    DivVis lp(vbd,H);
    
    Real a = 0.0;
    Real b[BL_SPACEDIM];
    b[0] = 1.0;
    b[1] = 1.0;
#if BL_SPACEDIM>2
    b[2] = 1.0;
#endif
    MultiFab  acoefs;
    int NcompA = (BL_SPACEDIM == 2  ?  2  :  1);
    acoefs.define(bs, NcompA, Nghost, Fab_allocate);
    acoefs.setVal(a);
    MultiFab bcoefs[BL_SPACEDIM];
    for (n=0; n<BL_SPACEDIM; ++n)
    {
	BoxArray bsC(bs);
	bcoefs[n].define(bsC.surroundingNodes(n), 1,
			 Nghost, Fab_allocate);
#if 1
	for(MFIter bmfi(bcoefs[n]); bmfi.isValid(); ++bmfi)
	{
	    FORT_MAKEMU(bcoefs[n][bmfi].dataPtr(),
			ARLIM(bcoefs[n][bmfi].loVect()),ARLIM(bcoefs[n][bmfi].hiVect()),H,n);
	}
#else
	bcoefs[n].setVal(b[n]);
#endif
    } // -->> over dimension
    lp.setCoefficients(acoefs, bcoefs);
#if 1
    lp.maxOrder(4);
#endif
    
    Nghost = 1;
    MultiFab tsoln(bs, Ncomp, Nghost, Fab_allocate); 
    tsoln.setVal(0.0);
#if 1
    tsoln.copy(fine);
#endif
#if 0
    // testing apply
    lp.apply(out,tsoln);
    Box subbox = out[0].box();
    Real n1 = out[0].norm(subbox,1,0,BL_SPACEDIM)*pow(H[0],BL_SPACEDIM);
    ParallelDescriptor::ReduceRealSum(n1);
    if (ParallelDescriptor::IOProcessor())
    {
	cout << "n1 output is "<<n1<<std::endl;
    }
    out.minus(rhs,0,BL_SPACEDIM,0);
    // special to single grid prob
    Real n2 = out[0].norm(subbox,1,0,BL_SPACEDIM)*pow(H[0],BL_SPACEDIM);
    ParallelDescriptor::ReduceRealSum(n2);
    if (ParallelDescriptor::IOProcessor())
    {
	cout << "n2 difference is "<<n2<<std::endl;
    }
#if 0
    subbox.grow(-1);
    Real n3 = out[0].norm(subbox,0,0,BL_SPACEDIM)*pow(H[0],BL_SPACEDIM);
    ParallelDescriptor::ReduceRealMax(n3);
    if (ParallelDescriptor::IOProcessor())
    {
	cout << "n3 difference is "<<n3<<std::endl;
    }
#endif
    
#endif
    
    const IntVect refRatio(D_DECL(2,2,2));
    const Real bgVal = 1.0;
    
#if 1
#ifndef NDEBUG
    // testing flux computation
    BoxArray xfluxbox(bs);
    xfluxbox.surroundingNodes(0);
    MultiFab xflux(xfluxbox,Ncomp,Nghost,Fab_allocate);
    xflux.setVal(1.e30);
    BoxArray yfluxbox(bs);
    yfluxbox.surroundingNodes(1);
    MultiFab yflux(yfluxbox,Ncomp,Nghost,Fab_allocate);
    yflux.setVal(1.e30);
#if BL_SPACEDIM>2
    BoxArray zfluxbox(bs);
    zfluxbox.surroundingNodes(2);
    MultiFab zflux(zfluxbox,Ncomp,Nghost,Fab_allocate);
    zflux.setVal(1.e30);
#endif
    lp.compFlux(xflux,
		yflux,
#if BL_SPACEDIM>2
		zflux,
#endif
		tsoln);
    
    // Write fluxes
    //writeMF(&xflux,"xflux.mfab");
    //writeMF(&yflux,"yflux.mfab");
#if BL_SPACEDIM>2
    //writeMF(&zflux,"zflux.mfab");
#endif
    
#endif
#endif
    
    Real tolerance = 1.0e-10; pp.query("tol", tolerance);
    Real tolerance_abs = 1.0e-10; pp.query("tol_abs", tolerance_abs);

#if 0
    cout << "Bndry Data object:" << std::endl;
    cout << lp.bndryData() << std::endl;
#endif
    
#if 0
    bool use_mg_pre = false;
    MCCGSolver cg(lp,use_mg_pre);
    cg.solve(soln,rhs,tolerance,tolerance_abs);
#else
    MCMultiGrid mg(lp);
    mg.solve(soln,rhs,tolerance,tolerance_abs);
#endif

#if 0
    cout << "MCLinOp object:" << std::endl;
    cout << lp << std::endl;
#endif
    
    VisMF::Write(soln,"soln");
    
#if 0
    // apply operator to soln to see if really satisfies eqn
    tsoln.copy(soln);
    lp.apply(out,tsoln);
    soln.copy(out);
    // Output "apply" results on soln
    VisMF::Write(soln,"apply");

    // Compute truncation
    for (MFIter smfi(soln); smfi.isValid(); ++smfi)
    {
	soln[smfi] -= fine[smfi];
    }
    for( int icomp=0; icomp < BL_SPACEDIM ; icomp++ )
    {
	Real solnMin = soln.min(icomp);
	Real solnMax = soln.max(icomp);
	ParallelDescriptor::ReduceRealMin(solnMin);
	ParallelDescriptor::ReduceRealMax(solnMax);
	if (ParallelDescriptor::IOProcessor())
	{
	    cout << icomp << "  "<<solnMin << " " << solnMax <<std::endl;
	}
    }
    // Output truncation
    VisMF::Write(soln,"trunc");
#endif

    int dumpLp=0; pp.query("dumpLp",dumpLp);
    bool write_lp = (dumpLp == 1 ? true : false);
    if (write_lp)
	std::cout << lp << std::endl;

    // Output trunc
    ParallelDescriptor::EndParallel();
}