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; } } }
void solve_with_Cpp(MultiFab& soln, MultiFab& gphi, Real a, Real b, MultiFab& alpha, PArray<MultiFab>& beta, MultiFab& rhs, const BoxArray& bs, const Geometry& geom) { BL_PROFILE("solve_with_Cpp()"); BndryData bd(bs, 1, geom); set_boundary(bd, rhs, 0); ABecLaplacian abec_operator(bd, dx); abec_operator.setScalars(a, b); abec_operator.setCoefficients(alpha, beta); MultiGrid mg(abec_operator); mg.setVerbose(verbose); mg.solve(soln, rhs, tolerance_rel, tolerance_abs); PArray<MultiFab> grad_phi(BL_SPACEDIM, PArrayManage); for (int n = 0; n < BL_SPACEDIM; ++n) grad_phi.set(n, new MultiFab(BoxArray(soln.boxArray()).surroundingNodes(n), 1, 0)); #if (BL_SPACEDIM == 2) abec_operator.compFlux(grad_phi[0],grad_phi[1],soln); #elif (BL_SPACEDIM == 3) abec_operator.compFlux(grad_phi[0],grad_phi[1],grad_phi[2],soln); #endif // Average edge-centered gradients to cell centers. BoxLib::average_face_to_cellcenter(gphi, grad_phi, geom); }
void solve_with_F90(MultiFab& soln, MultiFab& gphi, Real a, Real b, MultiFab& alpha, PArray<MultiFab>& beta, MultiFab& rhs, const BoxArray& bs, const Geometry& geom) { BL_PROFILE("solve_with_F90()"); FMultiGrid fmg(geom); int mg_bc[2*BL_SPACEDIM]; if (bc_type == Periodic) { // Define the type of boundary conditions to be periodic for ( int n = 0; n < BL_SPACEDIM; ++n ) { mg_bc[2*n + 0] = MGT_BC_PER; mg_bc[2*n + 1] = MGT_BC_PER; } } else if (bc_type == Neumann) { // Define the type of boundary conditions to be Neumann for ( int n = 0; n < BL_SPACEDIM; ++n ) { mg_bc[2*n + 0] = MGT_BC_NEU; mg_bc[2*n + 1] = MGT_BC_NEU; } } else if (bc_type == Dirichlet) { // Define the type of boundary conditions to be Dirichlet for ( int n = 0; n < BL_SPACEDIM; ++n ) { mg_bc[2*n + 0] = MGT_BC_DIR; mg_bc[2*n + 1] = MGT_BC_DIR; } } fmg.set_bc(mg_bc); fmg.set_maxorder(maxorder); fmg.set_scalars(a, b); fmg.set_coefficients(alpha, beta); int always_use_bnorm = 0; int need_grad_phi = 1; fmg.solve(soln, rhs, tolerance_rel, tolerance_abs, always_use_bnorm, need_grad_phi); PArray<MultiFab> grad_phi(BL_SPACEDIM, PArrayManage); for (int n = 0; n < BL_SPACEDIM; ++n) grad_phi.set(n, new MultiFab(BoxArray(soln.boxArray()).surroundingNodes(n), 1, 0)); fmg.get_fluxes(grad_phi); // Average edge-centered gradients to cell centers. BoxLib::average_face_to_cellcenter(gphi, grad_phi, geom); }
void ABec4::setCoefficients (const MultiFab &_a, const MultiFab &_b) { aCoefficients(_a); bCoefficients(_b); if (LO_Op) { int level = 0; const BoxArray& cba = boxArray(level); LO_Op->aCoefficients(_a); bool do_harm = true; for (int d=0; d<BL_SPACEDIM; ++d) { BoxArray eba = BoxArray(cba).surroundingNodes(d); MultiFab btmp(eba,1,0); lo_cc2ec(_b,btmp,0,0,1,d,do_harm); LO_Op->bCoefficients(btmp,d); } } }
void MCLinOp::prepareForLevel (int level) { if (level == 0) return; MCLinOp::prepareForLevel(level-1); if (h.size() > level) return; // // Assume from here down that this is a new level one coarser than existing // BL_ASSERT(h.size() == level); h.resize(level+1); int i; for (i = 0; i < BL_SPACEDIM; ++i) h[level][i] = h[level-1][i]*2.0; geomarray.resize(level+1); Box curdomain = Box( geomarray[level-1].Domain() ).coarsen(2); geomarray[level].define( curdomain ); // // Add a box to the new coarser level (assign removes old BoxArray) // gbox.resize(level+1); gbox[level] = BoxArray(gbox[level-1]).coarsen(2); // // Add the BndryRegister of relax values to the new coarser level. // BL_ASSERT(undrrelxr.size() == level); undrrelxr.resize(level+1); undrrelxr[level] = new BndryRegister(gbox[level], 1, 0, 0, numcomp); // // Add the BndryRegister to hold tagential derivatives to the new // coarser level. // BL_ASSERT(tangderiv.size() == level); tangderiv.resize(level+1); // // Figure out how many components. // const FabSet& samplefs = (*tangderiv[level-1])[Orientation(0,Orientation::low)]; tangderiv[level] = new BndryRegister(gbox[level],0,1,0,samplefs.nComp()); // // Add an Array of Array of maskvals to the new coarser level // For each orientation, build NULL masks, then use distributed allocation // Initial masks for coarse levels, ignore outside_domain possibility since // we always solve homogeneous equation on coarse levels. // BL_ASSERT(maskvals.size() == level); maskvals.resize(level+1); Array<IntVect> pshifts(27); std::vector< std::pair<int,Box> > isects; // // Use bgb's distribution map for masks. // We note that all orientations of the FabSets have the same distribution. // We'll use the low 0 side as the model. // for (FabSetIter bndryfsi(bgb[Orientation(0,Orientation::low)]); bndryfsi.isValid(); ++bndryfsi) { const int gn = bndryfsi.index(); MaskTuple& msk = maskvals[level][gn]; for (OrientationIter oitr; oitr; ++oitr) { const Orientation face = oitr(); Box bx_k = BoxLib::adjCell(gbox[level][gn], face, 1); // // Extend box in directions orthogonal to face normal. // for (int dir = 0; dir < BL_SPACEDIM; dir++) if (dir != face) bx_k.grow(dir,1); msk[face] = new Mask(bx_k, 1); Mask& curmask = *(msk[face]); curmask.setVal(BndryData::not_covered); gbox[level].intersections(bx_k,isects); for (int ii = 0, N = isects.size(); ii < N; ii++) if (isects[ii].first != gn) curmask.setVal(BndryData::covered, isects[ii].second, 0); // // Now take care of periodic wraparounds. // Geometry& curgeom = geomarray[level]; if (curgeom.isAnyPeriodic() && !curdomain.contains(bx_k)) { curgeom.periodicShift(curdomain, bx_k, pshifts); for (int iiv = 0, M = pshifts.size(); iiv < M; iiv++) { curmask.shift(pshifts[iiv]); gbox[level].intersections(curmask.box(),isects); for (int ii = 0, N = isects.size(); ii < N; ii++) curmask.setVal(BndryData::covered, isects[ii].second, 0); curmask.shift(-pshifts[iiv]); } } } } gbox[level].clear_hash_bin(); }
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 }
void solve_for_accel(PArray<MultiFab>& rhs, PArray<MultiFab>& phi, PArray<MultiFab>& grad_phi, const Array<Geometry>& geom, int base_level, int finest_level, Real offset) { Real tol = 1.e-10; Real abs_tol = 1.e-14; Array< PArray<MultiFab> > grad_phi_edge; grad_phi_edge.resize(rhs.size()); for (int lev = base_level; lev <= finest_level ; lev++) { grad_phi_edge[lev].resize(BL_SPACEDIM, PArrayManage); for (int n = 0; n < BL_SPACEDIM; ++n) grad_phi_edge[lev].set(n, new MultiFab(BoxArray(rhs[lev].boxArray()).surroundingNodes(n), 1, 1)); } Real strt = ParallelDescriptor::second(); // *************************************************** // Make sure the RHS sums to 0 if fully periodic // *************************************************** for (int lev = base_level; lev <= finest_level; lev++) { Real n0 = rhs[lev].norm0(); if (ParallelDescriptor::IOProcessor()) std::cout << "Max of rhs in solve_for_phi before correction at level " << lev << " " << n0 << std::endl; } for (int lev = base_level; lev <= finest_level; lev++) rhs[lev].plus(-offset, 0, 1, 0); for (int lev = base_level; lev <= finest_level; lev++) { Real n0 = rhs[lev].norm0(); if (ParallelDescriptor::IOProcessor()) std::cout << "Max of rhs in solve_for_phi after correction at level " << lev << " " << n0 << std::endl; } // *************************************************** // Solve for phi and return both phi and grad_phi_edge // *************************************************** #ifdef USEHPGMG solve_with_hpgmg(rhs,phi,grad_phi_edge,geom,base_level,finest_level,tol,abs_tol); #else solve_with_f90 (rhs,phi,grad_phi_edge,geom,base_level,finest_level,tol,abs_tol); #endif // Average edge-centered gradients to cell centers. for (int lev = base_level; lev <= finest_level; lev++) { BoxLib::average_face_to_cellcenter(grad_phi[lev], grad_phi_edge[lev], geom[lev]); geom[lev].FillPeriodicBoundary(grad_phi[lev],true); // wz: why only fill periodic boundary? } // VisMF::Write(grad_phi,"GradPhi"); { const int IOProc = ParallelDescriptor::IOProcessorNumber(); Real end = ParallelDescriptor::second() - strt; #if 0 #ifdef BL_LAZY Lazy::QueueReduction( [=] () mutable { #endif ParallelDescriptor::ReduceRealMax(end,IOProc); if (ParallelDescriptor::IOProcessor()) std::cout << "solve_for_phi() time = " << end << std::endl; #ifdef BL_LAZY }); #endif #endif } }
void solve_with_HPGMG(MultiFab& soln, MultiFab& gphi, Real a, Real b, MultiFab& alpha, PArray<MultiFab>& beta, MultiFab& beta_cc, MultiFab& rhs, const BoxArray& bs, const Geometry& geom, int n_cell) { BndryData bd(bs, 1, geom); set_boundary(bd, rhs, 0); ABecLaplacian abec_operator(bd, dx); abec_operator.setScalars(a, b); abec_operator.setCoefficients(alpha, beta); int minCoarseDim; if (domain_boundary_condition == BC_PERIODIC) { minCoarseDim = 2; // avoid problems with black box calculation of D^{-1} for poisson with periodic BC's on a 1^3 grid } else { minCoarseDim = 1; // assumes you can drop order on the boundaries } level_type level_h; mg_type MG_h; int numVectors = 12; int my_rank = 0, num_ranks = 1; #ifdef BL_USE_MPI MPI_Comm_size (MPI_COMM_WORLD, &num_ranks); MPI_Comm_rank (MPI_COMM_WORLD, &my_rank); #endif /* BL_USE_MPI */ const double h0 = dx[0]; // Create the geometric structure of the HPGMG grid using the RHS MultiFab as // a template. This doesn't copy any actual data. CreateHPGMGLevel(&level_h, rhs, n_cell, max_grid_size, my_rank, num_ranks, domain_boundary_condition, numVectors, h0); // Set up the coefficients for the linear operator L. SetupHPGMGCoefficients(a, b, alpha, beta_cc, &level_h); // Now that the HPGMG grid is built, populate it with RHS data. ConvertToHPGMGLevel(rhs, n_cell, max_grid_size, &level_h, VECTOR_F); #ifdef USE_HELMHOLTZ if (ParallelDescriptor::IOProcessor()) { std::cout << "Creating Helmholtz (a=" << a << ", b=" << b << ") test problem" << std::endl;; } #else if (ParallelDescriptor::IOProcessor()) { std::cout << "Creating Poisson (a=" << a << ", b=" << b << ") test problem" << std::endl;; } #endif /* USE_HELMHOLTZ */ if (level_h.boundary_condition.type == BC_PERIODIC) { double average_value_of_f = mean (&level_h, VECTOR_F); if (average_value_of_f != 0.0) { if (ParallelDescriptor::IOProcessor()) { std::cerr << "WARNING: Periodic boundary conditions, but f does not sum to zero... mean(f)=" << average_value_of_f << std::endl; } //shift_vector(&level_h,VECTOR_F,VECTOR_F,-average_value_of_f); } } //- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - rebuild_operator(&level_h,NULL,a,b); // i.e. calculate Dinv and lambda_max MGBuild(&MG_h,&level_h,a,b,minCoarseDim,ParallelDescriptor::Communicator()); // build the Multigrid Hierarchy //- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - if (ParallelDescriptor::IOProcessor()) std::cout << std::endl << std::endl << "===== STARTING SOLVE =====" << std::endl << std::flush; MGResetTimers (&MG_h); zero_vector (MG_h.levels[0], VECTOR_U); #ifdef USE_FCYCLES FMGSolve (&MG_h, 0, VECTOR_U, VECTOR_F, a, b, tolerance_abs, tolerance_rel); #else MGSolve (&MG_h, 0, VECTOR_U, VECTOR_F, a, b, tolerance_abs, tolerance_rel); #endif /* USE_FCYCLES */ MGPrintTiming (&MG_h, 0); // don't include the error check in the timing results //- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - if (ParallelDescriptor::IOProcessor()) std::cout << std::endl << std::endl << "===== Performing Richardson error analysis ==========================" << std::endl; // solve A^h u^h = f^h // solve A^2h u^2h = f^2h // solve A^4h u^4h = f^4h // error analysis... MGResetTimers(&MG_h); const double dtol = tolerance_abs; const double rtol = tolerance_rel; int l;for(l=0;l<3;l++){ if(l>0)restriction(MG_h.levels[l],VECTOR_F,MG_h.levels[l-1],VECTOR_F,RESTRICT_CELL); zero_vector(MG_h.levels[l],VECTOR_U); #ifdef USE_FCYCLES FMGSolve(&MG_h,l,VECTOR_U,VECTOR_F,a,b,dtol,rtol); #else MGSolve(&MG_h,l,VECTOR_U,VECTOR_F,a,b,dtol,rtol); #endif } richardson_error(&MG_h,0,VECTOR_U); // Now convert solution from HPGMG back to rhs MultiFab. ConvertFromHPGMGLevel(soln, &level_h, VECTOR_U); const double norm_from_HPGMG = norm(&level_h, VECTOR_U); const double mean_from_HPGMG = mean(&level_h, VECTOR_U); const Real norm0 = soln.norm0(); const Real norm2 = soln.norm2(); if (ParallelDescriptor::IOProcessor()) { std::cout << "mean from HPGMG: " << mean_from_HPGMG << std::endl; std::cout << "norm from HPGMG: " << norm_from_HPGMG << std::endl; std::cout << "norm0 of RHS copied to MF: " << norm0 << std::endl; std::cout << "norm2 of RHS copied to MF: " << norm2 << std::endl; } // Write the MF to disk for comparison with the in-house solver if (plot_soln) { writePlotFile("SOLN-HPGMG", soln, geom); } //- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - MGDestroy(&MG_h); destroy_level(&level_h); //- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - PArray<MultiFab> grad_phi(BL_SPACEDIM, PArrayManage); for (int n = 0; n < BL_SPACEDIM; ++n) grad_phi.set(n, new MultiFab(BoxArray(soln.boxArray()).surroundingNodes(n), 1, 0)); #if (BL_SPACEDIM == 2) abec_operator.compFlux(grad_phi[0],grad_phi[1],soln); #elif (BL_SPACEDIM == 3) abec_operator.compFlux(grad_phi[0],grad_phi[1],grad_phi[2],soln); #endif // Average edge-centered gradients to cell centers. BoxLib::average_face_to_cellcenter(gphi, grad_phi, geom); }
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(); }