forceinline int ComplementView<View>::glbMin(void) const { LubRanges<View> ub(x); RangesCompl<LubRanges<View> > ubc(ub); if (ubc()) { return ubc.min(); } else { return BndSet::MIN_OF_EMPTY; } }
forceinline int ComplementView<View>::glbMax(void) const { LubRanges<View> ub(x); RangesCompl<LubRanges<View> > ubc(ub); if (ubc()) { while (ubc()) ++ubc; return ubc.max(); } else { return BndSet::MAX_OF_EMPTY; } }
/** * Probably should be pushed back into ContField? */ void IncNavierStokes::SetRadiationBoundaryForcing(int fieldid) { int i,n; Array<OneD, const SpatialDomains::BoundaryConditionShPtr > BndConds; Array<OneD, MultiRegions::ExpListSharedPtr> BndExp; BndConds = m_fields[fieldid]->GetBndConditions(); BndExp = m_fields[fieldid]->GetBndCondExpansions(); StdRegions::StdExpansionSharedPtr elmt; StdRegions::StdExpansionSharedPtr Bc; int cnt; int elmtid,nq,offset, boundary; Array<OneD, NekDouble> Bvals, U; int cnt1 = 0; for(cnt = n = 0; n < BndConds.num_elements(); ++n) { std::string type = BndConds[n]->GetUserDefined(); if((BndConds[n]->GetBoundaryConditionType() == SpatialDomains::eRobin)&&(boost::iequals(type,"Radiation"))) { for(i = 0; i < BndExp[n]->GetExpSize(); ++i,cnt++) { elmtid = m_fieldsBCToElmtID[fieldid][cnt]; elmt = m_fields[fieldid]->GetExp(elmtid); offset = m_fields[fieldid]->GetPhys_Offset(elmtid); U = m_fields[fieldid]->UpdatePhys() + offset; Bc = BndExp[n]->GetExp(i); boundary = m_fieldsBCToTraceID[fieldid][cnt]; // Get edge values and put into ubc nq = Bc->GetTotPoints(); Array<OneD, NekDouble> ubc(nq); elmt->GetTracePhysVals(boundary,Bc,U,ubc); Vmath::Vmul(nq,&m_fieldsRadiationFactor[fieldid][cnt1 + BndExp[n]->GetPhys_Offset(i)],1,&ubc[0],1,&ubc[0],1); Bvals = BndExp[n]->UpdateCoeffs()+BndExp[n]->GetCoeff_Offset(i); Bc->IProductWRTBase(ubc,Bvals); } cnt1 += BndExp[n]->GetTotPoints(); } else { cnt += BndExp[n]->GetExpSize(); } } }
//--------------------------------------------------------- void TestPoissonIPDG3D::Run() //--------------------------------------------------------- { umLOG(1, "TestPoissonIPDG3D::Run()\n"); double wt1=timer.read(); // Initialize solver and construct grid and metric StartUp3D(); #if (0) //------------------------------------- // check the mesh //------------------------------------- Output_Mesh(); umLOG(1, "\n*** Exiting after writing mesh\n\n"); return; #endif // sparse operators CSd A("OP"), M("MM"); // build 3D IPDG Poisson matrix (assuming all Dirichlet) PoissonIPDG3D(A, M); if (0) { // NBN: experiment with diagonal strength int Nz=0; for (int j=0; j<A.n; ++j) { for (int p = A.P[j]; p<A.P[j+1]; ++p) { if (A.I[p] == j) { A.X[p] += A.X[p]; // augment A(i,j) } } } } if (0) { umLOG(1, "Dumping file Afull.dat for Matlab\n"); umLOG(1, "A(%d,%d), nnz = %d\n", A.m,A.n, A.P[A.n]); FILE* fp=fopen("Afull.dat", "w"); A.write_ML(fp); fclose(fp); umLOG(1, "exiting.\n"); return; } // iterative solver CS_PCG it_sol; try { // Note: operator A is symmetric, with only its // lower triangule stored. We use this symmetry // to accelerate operator A*x int flag=sp_SYMMETRIC; flag|=sp_LOWER; flag|=sp_TRIANGULAR; A.set_shape(flag); // drop tolerance for cholinc double droptol=1e-4; // Note: ownership of A is transfered to solver object it_sol.cholinc(A, droptol); } catch(...) { umLOG(1, "\nCaught exception from symbolic chol.\n"); } DVec exact("exact"), f("f"), rhs("rhs"), u("u"); double t1=0.0,t2=0.0; //------------------------------------------------------- // set up boundary condition //------------------------------------------------------- DVec xbc = Fx(mapB), ybc = Fy(mapB), zbc = Fz(mapB); DVec ubc(Nfp*Nfaces*K); ubc(mapB) = sin(pi*xbc) * sin(pi*ybc) * sin(pi*zbc); //------------------------------------------------------- // form right hand side contribution from boundary condition //------------------------------------------------------- DMat Abc = PoissonIPDGbc3D(ubc); //------------------------------------------------------- // evaluate forcing function //------------------------------------------------------- exact = sin(pi*x).dm(sin(pi*y).dm(sin(pi*z))); f = (-3.0*SQ(pi)) * exact; //------------------------------------------------------- // set up right hand side for variational Poisson equation //------------------------------------------------------- rhs = M*(-f) + (DVec&)(Abc); //------------------------------------------------------- // solve using pcg iterative solver //------------------------------------------------------- t1 = timer.read(); u = it_sol.solve(rhs, 1e-9, 30); t2 = timer.read(); //------------------------------------------------------- // compute nodal error //------------------------------------------------------- r = (u-exact); m_maxAbsError = r.max_val_abs(); umLOG(1, " solve -- done. (%0.4lf sec)\n\n", t2-t1); umLOG(1, " max error = %g\n\n", m_maxAbsError); #if (0) //####################################################### // plot solution and error figure(2); subplot(1,3,2); PlotContour3D(2*N, u, linspace(-1, 1, 10)); title('numerical solution'); subplot(1,3,3); PlotContour3D(2*N, log10(eps+abs(u-exact)), linspace(-4, 0, 10)); title('numerical error'); //####################################################### #endif double wt2=timer.read(); umLOG(1, "TestPoissonIPDG3D::Run() complete\n"); umLOG(1, "total time = %0.4lf sec\n\n", wt2-wt1); }
//--------------------------------------------------------- DVec& NDG2D::PoissonIPDGbc2D (DVec& ubc, //[in] DVec& qbc //[in] ) //--------------------------------------------------------- { // function [OP] = PoissonIPDGbc2D() // Purpose: Set up the discrete Poisson matrix directly // using LDG. The operator is set up in the weak form // build DG derivative matrices int max_OP = (K*Np*Np*(1+Nfaces)); // initialize parameters DVec faceR("faceR"), faceS("faceS"); DMat V1D("V1D"), Dx("Dx"),Dy("Dy"), Dn1("Dn1"), mmE_Fm1("mmE(:,Fm1)"); IVec Fm("Fm"), Fm1("Fm1"), fidM("fidM"); double lnx=0.0,lny=0.0,lsJ=0.0,hinv=0.0,gtau=0.0; int i=0,k1=0,f1=0,id=0; IVec i1_Nfp = Range(1,Nfp); double N1N1 = double((N+1)*(N+1)); // build local face matrices DMat massEdge[4]; // = zeros(Np,Np,Nfaces); for (i=1; i<=Nfaces; ++i) { massEdge[i].resize(Np,Np); } // face mass matrix 1 Fm = Fmask(All,1); faceR = r(Fm); V1D = Vandermonde1D(N, faceR); massEdge[1](Fm,Fm) = inv(V1D*trans(V1D)); // face mass matrix 2 Fm = Fmask(All,2); faceR = r(Fm); V1D = Vandermonde1D(N, faceR); massEdge[2](Fm,Fm) = inv(V1D*trans(V1D)); // face mass matrix 3 Fm = Fmask(All,3); faceS = s(Fm); V1D = Vandermonde1D(N, faceS); massEdge[3](Fm,Fm) = inv(V1D*trans(V1D)); // build DG right hand side DVec* pBC = new DVec(Np*K, "bc", OBJ_temp); DVec& bc = (*pBC); // reference, for syntax //////////////////////////////////////////////////////////////// umMSG(1, "\n ==> {OP} assembly [bc]: "); for (k1=1; k1<=K; ++k1) { if (! (k1%100)) { umMSG(1, "%d, ",k1); } // rows1 = outer(Range((k1-1)*Np+1,k1*Np), Ones(NGauss)); // Build element-to-element parts of operator for (f1=1; f1<=Nfaces; ++f1) { if (BCType(k1,f1)) { ////////////////////////added by Kevin /////////////////////////////// Fm1 = Fmask(All,f1); fidM = (k1-1)*Nfp*Nfaces + (f1-1)*Nfp + i1_Nfp; id = 1+(f1-1)*Nfp + (k1-1)*Nfp*Nfaces; lnx = nx(id); lny = ny(id); lsJ = sJ(id); hinv = Fscale(id); Dx = rx(1,k1)*Dr + sx(1,k1)*Ds; Dy = ry(1,k1)*Dr + sy(1,k1)*Ds; Dn1 = lnx*Dx + lny*Dy; //mmE = lsJ*massEdge(:,:,f1); //bc(All,k1) += (gtau*mmE(All,Fm1) - Dn1'*mmE(All,Fm1))*ubc(fidM); mmE_Fm1 = massEdge[f1](All,Fm1); mmE_Fm1 *= lsJ; gtau = 10*N1N1*hinv; // set penalty scaling //bc(All,k1) += (gtau*mmE_Fm1 - trans(Dn1)*mmE_Fm1) * ubc(fidM); switch(BCType(k1,f1)){ case BC_Dirichlet: bc(Np*(k1-1)+Range(1,Np)) += (gtau*mmE_Fm1 - trans(Dn1)*mmE_Fm1)*ubc(fidM); break; case BC_Neuman: bc(Np*(k1-1)+Range(1,Np)) += mmE_Fm1*qbc(fidM); break; default: std::cout<<"warning: boundary condition is incorrect"<<std::endl; } } } } return bc; }