void ConvectionDiffusion::DiffusionFluxes(abstract_variable & param, Face fKL, variable & fluxT, variable & corrK, variable & corrL, variable & corrK_cK, variable & corrL_cL) const { if( !BuildFlux(fKL) ) return; get_variable<variable> Func(param); //diffusion part if( tag_K.isValid() ) { Cell cK = fKL->BackCell(); Cell cL = fKL->FrontCell(); //diffusion data real coefT = fKL.Real (tag_DIFF_CT); real rhsT = fKL.Real (tag_DIFF_RT); //corrector in back cell real_array coefsK = fKL->RealArray (tag_DIFF_CB); ref_array elemsK = fKL->ReferenceArray(tag_DIFF_EB); real rhsK = fKL->Real (tag_DIFF_RB); //corrector in front cell real_array coefsL = fKL->RealArray (tag_DIFF_CF); ref_array elemsL = fKL->ReferenceArray(tag_DIFF_EF); real rhsL = fKL->Real (tag_DIFF_RF); //flag bulk flag = fKL->Bulk (tag_DIFF_F); //two-point part fluxT -= rhsT; if( flag > 1 ) fluxT -= coefT*(Func(cL) - Func(cK)); else fluxT += coefT*Func(cK); if( flag > 1 ) //full nonlinear part { if( perform_correction_diffusion ) { corrK -= rhsK; corrL -= rhsL; for(enumerator k = 0; k < elemsK.size(); ++k) if( elemsK[k].isValid() ) corrK -= Func(elemsK[k])*coefsK[k]; for(enumerator k = 0; k < elemsL.size(); ++k) if( elemsL[k].isValid() ) corrL -= Func(elemsL[k])*coefsL[k]; //if( scheme_type == NTPFA || scheme_type == NTPFA_PICARD ) { if( elemsK.back().isValid() ) corrK_cK -= Func(elemsK.back())*coefsK.back(); if( elemsL.back().isValid() ) corrL_cL -= Func(elemsL.back())*coefsL.back(); } } } else //one-sided corrector { if( flag == 1 && perform_correction_diffusion ) // only fluxK is present { corrK -= rhsK; for(enumerator k = 0; k < elemsK.size(); ++k) if( elemsK[k].isValid() ) corrK -= Func(elemsK[k])*coefsK[k]; } //no corrections otherwise } //flag } //diffusion part }
int main(int argc,char ** argv) { Solver::Initialize(&argc,&argv,""); // Initialize the solver and MPI activity #if defined(USE_PARTITIONER) Partitioner::Initialize(&argc,&argv); // Initialize the partitioner activity #endif if( argc > 1 ) { TagReal phi; TagReal tag_F; TagRealArray tag_K; TagRealArray tag_BC; TagReal phi_ref; Mesh * m = new Mesh(); // Create an empty mesh double ttt = Timer(); bool repartition = false; m->SetCommunicator(INMOST_MPI_COMM_WORLD); // Set the MPI communicator for the mesh if( m->GetProcessorRank() == 0 ) // If the current process is the master one std::cout << argv[0] << std::endl; if( m->isParallelFileFormat(argv[1]) ) { m->Load(argv[1]); // Load mesh from the parallel file format repartition = true; } else { if( m->GetProcessorRank() == 0 ) m->Load(argv[1]); // Load mesh from the serial file format } BARRIER; if( m->GetProcessorRank() == 0 ) std::cout << "Processors: " << m->GetProcessorsNumber() << std::endl; if( m->GetProcessorRank() == 0 ) std::cout << "Load(MPI_File): " << Timer()-ttt << std::endl; //~ double ttt2 = Timer(); //~ Mesh t; //~ t.SetCommunicator(INMOST_MPI_COMM_WORLD); //~ t.SetParallelFileStrategy(0); //~ t.Load(argv[1]); //~ BARRIER //~ if( m->GetProcessorRank() == 0 ) std::cout << "Load(MPI_Scatter): " << Timer()-ttt2 << std::endl; #if defined(USE_PARTITIONER) if (m->GetProcessorsNumber() > 1) { // currently only non-distributed meshes are supported by Inner_RCM partitioner ttt = Timer(); Partitioner * p = new Partitioner(m); p->SetMethod(Partitioner::INNER_KMEANS,Partitioner::Partition); // Specify the partitioner p->Evaluate(); // Compute the partitioner and store new processor ID in the mesh delete p; BARRIER; if( m->GetProcessorRank() == 0 ) std::cout << "Evaluate: " << Timer()-ttt << std::endl; ttt = Timer(); m->Redistribute(); // Redistribute the mesh data m->ReorderEmpty(CELL|FACE|EDGE|NODE); // Clean the data after reordring BARRIER; if( m->GetProcessorRank() == 0 ) std::cout << "Redistribute: " << Timer()-ttt << std::endl; } #endif ttt = Timer(); phi = m->CreateTag("Solution",DATA_REAL,CELL,NONE,1); // Create a new tag for the solution phi bool makerefsol = true; if( m->HaveTag("PERM" ) ) { tag_K = m->GetTag("PERM"); makerefsol = false; std::cout << "Permeability from grid" << std::endl; } else { std::cout << "Set perm" << std::endl; tag_K = m->CreateTag("PERM",DATA_REAL,CELL,NONE,1); // Create a new tag for K tensor for( Mesh::iteratorCell cell = m->BeginCell(); cell != m->EndCell(); ++cell ) // Loop over mesh cells tag_K[*cell][0] = 1.0; // Store the tensor K value into the tag } if( m->HaveTag("BOUNDARY_CONDITION") ) { tag_BC = m->GetTag("BOUNDARY_CONDITION"); makerefsol = false; std::cout << "Boundary conditions from grid" << std::endl; } else { std::cout << "Set boundary conditions" << std::endl; double x[3]; tag_BC = m->CreateTag("BOUNDARY_CONDITION",DATA_REAL,FACE,FACE,3); for( Mesh::iteratorFace face = m->BeginFace(); face != m->EndFace(); ++face ) if( face->Boundary() && !(face->GetStatus() == Element::Ghost) ) { face->Centroid(x); tag_BC[*face][0] = 1; //dirichlet tag_BC[*face][1] = 0; //neumann tag_BC[*face][2] = func(x,0);//face->Mean(func, 0); } } if( m->HaveTag("FORCE") ) { tag_F = m->GetTag("FORCE"); makerefsol = false; std::cout << "Force from grid" << std::endl; } else if( makerefsol ) { std::cout << "Set rhs" << std::endl; tag_F = m->CreateTag("FORCE",DATA_REAL,CELL,NONE,1); // Create a new tag for external force double x[3]; for( Mesh::iteratorCell cell = m->BeginCell(); cell != m->EndCell(); ++cell ) // Loop over mesh cells { cell->Centroid(x); tag_F[*cell] = -func_rhs(x,1); //tag_F[*cell] = -cell->Mean(func_rhs,1); } } if(m->HaveTag("REFERENCE_SOLUTION") ) phi_ref = m->GetTag("REFERENCE_SOLUTION"); else if( makerefsol ) { phi_ref = m->CreateTag("REFRENCE_SOLUTION",DATA_REAL,CELL,NONE,1); double x[3]; for( Mesh::iteratorCell cell = m->BeginCell(); cell != m->EndCell(); ++cell ) { cell->Centroid(x); phi_ref[*cell] = func(x,0);//cell->Mean(func, 0); } } ttt = Timer(); m->ExchangeGhost(1,FACE); m->ExchangeData(tag_K,CELL,0); // Exchange the tag_K data over processors BARRIER; if( m->GetProcessorRank() == 0 ) std::cout << "Exchange ghost: " << Timer()-ttt << std::endl; ttt = Timer(); Solver S("inner_ilu2"); // Specify the linear solver to ASM+ILU2+BiCGStab one S.SetParameter("absolute_tolerance", "1e-8"); S.SetParameter("schwartz_overlap", "2"); Residual R; // Residual vector Sparse::LockService Locks; Sparse::Vector Update; // Declare the solution and the right-hand side vectors { Mesh::GeomParam table; table[CENTROID] = CELL | FACE; table[NORMAL] = FACE; table[ORIENTATION] = FACE; table[MEASURE] = CELL | FACE; table[BARYCENTER] = CELL | FACE; m->PrepareGeometricData(table); } BARRIER if( m->GetProcessorRank() == 0 ) std::cout << "Prepare geometric data: " << Timer()-ttt << std::endl; { Automatizator aut; Automatizator::MakeCurrent(&aut); INMOST_DATA_ENUM_TYPE iphi = aut.RegisterTag(phi,CELL); aut.EnumerateEntries(); // Set the indeces intervals for the matrix and vectors R.SetInterval(aut.GetFirstIndex(),aut.GetLastIndex()); R.InitLocks(); Update.SetInterval(aut.GetFirstIndex(),aut.GetLastIndex()); dynamic_variable Phi(aut,iphi); // Solve \nabla \cdot \nabla phi = f equation //for( Mesh::iteratorFace face = m->BeginFace(); face != m->EndFace(); ++face ) #if defined(USE_OMP) #pragma omp parallel #endif { variable flux; //should be more efficient to define here to avoid multiple memory allocations if storage for variations should be expanded rMatrix x1(3,1), x2(3,1), xf(3,1), n(3,1); double d1, d2, k1, k2, area, T, a, b, c; #if defined(USE_OMP) #pragma omp for #endif for(Storage::integer iface = 0; iface < m->FaceLastLocalID(); ++iface ) if( m->isValidFace(iface) ) { Face face = Face(m,ComposeFaceHandle(iface)); Element::Status s1,s2; Cell r1 = face->BackCell(); Cell r2 = face->FrontCell(); if( ((!r1->isValid() || (s1 = r1->GetStatus()) == Element::Ghost)?0:1) + ((!r2->isValid() || (s2 = r2->GetStatus()) == Element::Ghost)?0:1) == 0) continue; area = face->Area(); // Get the face area face->UnitNormal(n.data()); // Get the face normal face->Centroid(xf.data()); // Get the barycenter of the face r1->Centroid(x1.data()); // Get the barycenter of the cell k1 = n.DotProduct(rMatrix::FromTensor(tag_K[r1].data(), tag_K[r1].size(),3)*n); d1 = fabs(n.DotProduct(xf-x1)); if( !r2->isValid() ) // boundary condition { // bnd_pnt is a projection of the cell center to the face // a*pb + bT(pb-p1) = c // F = T(pb-p1) // pb = (c + bTp1)/(a+bT) // F = T/(a+bT)(c - ap1) T = k1/d1; a = 0; b = 1; c = 0; if( tag_BC.isValid() && face.HaveData(tag_BC) ) { a = tag_BC[face][0]; b = tag_BC[face][1]; c = tag_BC[face][2]; //std::cout << "a " << a << " b " << b << " c " << c << std::endl; } R.Lock(Phi.Index(r1)); R[Phi.Index(r1)] -= T/(a + b*T) * area * (c - a*Phi(r1)); R.UnLock(Phi.Index(r1)); } else { r2->Centroid(x2.data()); k2 = n.DotProduct(rMatrix::FromTensor(tag_K[r2].data(), tag_K[r2].size(),3)*n); d2 = fabs(n.DotProduct(x2-xf)); T = 1.0/(d1/k1 + d2/k2); flux = T* area * (Phi(r2) - Phi(r1)); if( s1 != Element::Ghost ) { R.Lock(Phi.Index(r1)); R[Phi.Index(r1)] -= flux; R.UnLock(Phi.Index(r1)); } if( s2 != Element::Ghost ) { R.Lock(Phi.Index(r2)); R[Phi.Index(r2)] += flux; R.UnLock(Phi.Index(r2)); } } } } if( tag_F.isValid() ) { #if defined(USE_OMP) #pragma omp parallel for #endif for( Storage::integer icell = 0; icell < m->CellLastLocalID(); ++icell ) if( m->isValidCell(icell) ) { Cell cell = Cell(m,ComposeCellHandle(icell)); if( cell->GetStatus() != Element::Ghost ) R[Phi.Index(cell)] -= tag_F[cell] * cell->Volume(); } } BARRIER; if( m->GetProcessorRank() == 0 ) std::cout << "Matrix assemble: " << Timer()-ttt << std::endl; //m->RemoveGeometricData(table); // Clean the computed geometric data if( argc > 3 ) // Save the matrix and RHS if required { ttt = Timer(); R.GetJacobian().Save(std::string(argv[2])); // "A.mtx" R.GetResidual().Save(std::string(argv[3])); // "b.rhs" BARRIER; if( m->GetProcessorRank() == 0 ) std::cout << "Save matrix \"" << argv[2] << "\" and RHS \"" << argv[3] << "\": " << Timer()-ttt << std::endl; } ttt = Timer(); S.SetMatrix(R.GetJacobian()); // Compute the preconditioner for the original matrix S.Solve(R.GetResidual(),Update); // Solve the linear system with the previously computted preconditioner BARRIER; if( m->GetProcessorRank() == 0 ) { std::cout << S.Residual() << " " << S.Iterations() << " " << S.ReturnReason() << std::endl; std::cout << "Solve system: " << Timer()-ttt << std::endl; } ttt = Timer(); if( phi_ref.isValid() ) { Tag error = m->CreateTag("error",DATA_REAL,CELL,NONE,1); double err_C = 0.0, err_L2 = 0.0, vol = 0.0; #if defined(USE_OMP) #pragma omp parallel #endif { double local_err_C = 0; #if defined(USE_OMP) #pragma omp for reduction(+:err_L2) reduction(+:vol) #endif for( Storage::integer icell = 0; icell < m->CellLastLocalID(); ++icell ) if( m->isValidCell(icell) ) { Cell cell = Cell(m,ComposeCellHandle(icell)); if( cell->GetStatus() != Element::Ghost ) { double old = phi[cell]; double exact = phi_ref[cell]; double res = Update[Phi.Index(cell)]; double sol = old-res; double err = fabs (sol - exact); if (err > local_err_C) local_err_C = err; err_L2 += err * err * cell->Volume(); vol += cell->Volume(); cell->Real(error) = err; phi[cell] = sol; } } #if defined(USE_OMP) #pragma omp critical #endif { if( local_err_C > err_C ) err_C = local_err_C; } } err_C = m->AggregateMax(err_C); // Compute the maximal C norm for the error err_L2 = sqrt(m->Integrate(err_L2)/m->Integrate(vol)); // Compute the global L2 norm for the error if( m->GetProcessorRank() == 0 ) std::cout << "err_C = " << err_C << std::endl; if( m->GetProcessorRank() == 0 ) std::cout << "err_L2 = " << err_L2 << std::endl; } } BARRIER; if( m->GetProcessorRank() == 0 ) std::cout << "Compute true residual: " << Timer()-ttt << std::endl; ttt = Timer(); m->ExchangeData(phi,CELL,0); // Data exchange over processors BARRIER; if( m->GetProcessorRank() == 0 ) std::cout << "Exchange phi: " << Timer()-ttt << std::endl; std::string filename = "result"; if( m->GetProcessorsNumber() == 1 ) filename += ".vtk"; else filename += ".pvtk"; ttt = Timer(); m->Save(filename); m->Save("result.pmf"); BARRIER; if( m->GetProcessorRank() == 0 ) std::cout << "Save \"" << filename << "\": " << Timer()-ttt << std::endl; delete m; }
ConvectionDiffusion::ConvectionDiffusion(Mesh * _m, Tag _tag_U, Tag _tag_K, Tag _tag_BC, MarkerType _boundary_face, bool correct_convection, bool correct_diffusion) : m(_m), tag_U(_tag_U), tag_K(_tag_K), tag_BC(_tag_BC), boundary_face(_boundary_face) { perform_correction_diffusion = correct_convection; perform_correction_convection = correct_diffusion; if( tag_U.isValid() ) { //4 entries - triplet for upwind corrector, including the cell //back cell corrector tag_CONV_CB = m->CreateTag("CONV_BACK_COEFS" ,DATA_REAL ,FACE,NONE,4); tag_CONV_EB = m->CreateTag("CONV_BACK_ELEMS" ,DATA_REFERENCE,FACE,NONE,4); tag_CONV_RB = m->CreateTag("CONV_BACK_RHS" ,DATA_REAL ,FACE,NONE,1); //front cell corrector tag_CONV_CF = m->CreateTag("CONV_FRONT_COEFS" ,DATA_REAL ,FACE,NONE,4); tag_CONV_EF = m->CreateTag("CONV_FRONT_ELEMS" ,DATA_REFERENCE,FACE,NONE,4); tag_CONV_RF = m->CreateTag("CONV_FRONT_RHS" ,DATA_REAL ,FACE,NONE,1); //upstream tag_CONV_CU = m->CreateTag("CONV_UPSTREAM_COEF",DATA_REAL ,FACE,NONE,1); tag_CONV_EU = m->CreateTag("CONV_UPSTREAM_ELEM",DATA_REFERENCE,FACE,NONE,1); tag_CONV_RU = m->CreateTag("CONV_UPSTREAM_RHS" ,DATA_REAL ,FACE,NONE,1); tag_CONV_VU = m->CreateTag("CONV_UPSTREAM_CORR",DATA_REAL ,FACE,NONE,6); //flag tag_CONV_F = m->CreateTag("CONV_FLAG" ,DATA_BULK ,FACE,NONE,1); } if( tag_K.isValid() ) { //4 entries - triplet for transversal corrector, including the cell //back cell corrector tag_DIFF_CB = m->CreateTag("DIFF_BACK_COEFS" ,DATA_REAL ,FACE,NONE,4); tag_DIFF_EB = m->CreateTag("DIFF_BACK_ELEMS" ,DATA_REFERENCE,FACE,NONE,4); tag_DIFF_RB = m->CreateTag("DIFF_BACK_RHS" ,DATA_REAL ,FACE,NONE,1); //front cell corrector tag_DIFF_CF = m->CreateTag("DIFF_FRONT_COEFS" ,DATA_REAL ,FACE,NONE,4); tag_DIFF_EF = m->CreateTag("DIFF_FRONT_ELEMS" ,DATA_REFERENCE,FACE,NONE,4); tag_DIFF_RF = m->CreateTag("DIFF_FRONT_RHS" ,DATA_REAL ,FACE,NONE,1); //two-point transmissibility tag_DIFF_CT = m->CreateTag("DIFF_TP_COEF" ,DATA_REAL ,FACE,NONE,1); tag_DIFF_RT = m->CreateTag("DIFF_TP_RHS" ,DATA_REAL ,FACE,NONE,1); tag_DIFF_VT = m->CreateTag("DIFF_TP_CORR" ,DATA_REAL ,FACE,NONE,6); //flag tag_DIFF_F = m->CreateTag("DIFF_FLAG" ,DATA_BULK ,FACE,NONE,1); } build_flux = 0; if( m->GetProcessorsNumber() > 1 ) { build_flux = m->CreateMarker(); #if defined(USE_OMP) #pragma omp for #endif for(integer q = 0; q < m->FaceLastLocalID(); ++q ) if( m->isValidFace(q) ) { Face fKL = m->FaceByLocalID(q); Element::Status stat = fKL->BackCell()->GetStatus(); if( fKL->FrontCell().isValid() ) stat |= fKL->FrontCell()->GetStatus(); if( stat & (Element::Shared | Element::Owned) ) fKL->SetMarker(build_flux); } } initialization_failure = false; //there is a failure during intialization //Precompute two-point part, transversal correction direction, // upstream and upstream correction directions and flags #if defined(USE_OMP) #pragma omp parallel #endif { //private memory rMatrix xK(1,3), //center of cell K yK(1,3), //projection of xK onto interface xL(1,3), //center of cell L yL(1,3), //projection of xL onto interface xKL(1,3), //center of face common to K and L nKL(1,3), //normal vector to face KKn(1,3), //co-normal at cell K KLn(1,3), //co-normal at cell L v(1,3), //vector for correction uK(1,3), //vector for upwind correction in back cell uL(1,3), //vector for upwind correction in front cell r(1,3), //path to reach upstream concentration from downstream concentration r0(1,3), KK(3,3), //tenosr at cell K KL(3,3), //tensor at cell L KD(3,3), //difference of tensors gammaK(1,3), //transversal part of co-normal at cell K gammaL(1,3), //transversal part of co-normal at cell L iT(3,3), //heterogeneous interpolation tensor iC(1,3) //heterogeneous interpolation correction ; real A, //area of the face U, //normal component of the velocity C, //coefficient for upstream cell T, //two-point transmissibility R, //right hand side dK, //distance from center to interface at cell K dL, //distance from center to interface at cell L lambdaK, //projection of co-normal onto normal at cell K lambdaL //projection of co-normal onto normal at cell L ; const real eps = degenerate_diffusion_regularization; Cell cK, cL, cU; Face fKL; bulk flag_DIFF, flag_CONV; KK.Zero(); KL.Zero(); KD.Zero(); U = 0.0; #if defined(USE_OMP) #pragma omp for #endif for(integer q = 0; q < m->FaceLastLocalID(); ++q ) if( m->isValidFace(q) ) { fKL = m->FaceByLocalID(q); if( !BuildFlux(fKL) ) continue; fKL.Centroid(xKL.data()); fKL.UnitNormal(nKL.data()); A = fKL.Area(); if( tag_U.isValid() ) U = fKL.Real(tag_U); cK = fKL.BackCell(); cL = fKL.FrontCell(); assert(cK.isValid()); cK.Centroid(xK.data()); dK = nKL.DotProduct(xKL-xK); yK = xK + dK*nKL; if( tag_K.isValid() ) KK = rMatrix::FromTensor(cK.RealArray(tag_K).data(),cK.RealArray(tag_K).size());//.Transpose(); KKn = nKL*KK; lambdaK = nKL.DotProduct(KKn); gammaK = KKn - lambdaK*nKL; //Diffusion part uK.Zero(); uL.Zero(); if( cL.isValid() ) //internal, both cells are present { cL.Centroid(xL.data()); dL = nKL.DotProduct(xL-xKL); yL = xL - dL*nKL; if( tag_K.isValid() ) KL = rMatrix::FromTensor(cL.RealArray(tag_K).data(),cL.RealArray(tag_K).size());//.Transpose(); KLn = nKL*KL; lambdaL = nKL.DotProduct(KLn); gammaL = KLn - lambdaL*nKL; //don't forget the area! R = 0; if( split_diffusion ) { uK = A * (lambdaK*lambdaL*(yK-yL) + lambdaL*dK*gammaK + lambdaK*dL*gammaL)/(dK*lambdaL+dL*lambdaK + 1.0e-100); uL = uK; T = A * lambdaK*lambdaL/(dK*lambdaL+dL*lambdaK + 1.0e-100); } else //un-splitted diffusion { uK = A * KKn; uL = A * KLn; T = 0; } //internal if( uK.FrobeniusNorm() > 0 || uL.FrobeniusNorm() > 0 ) flag_DIFF = 3; else flag_DIFF = 2; } else if( fKL.GetMarker(boundary_face) ) //boundary interface { //don't forget the area! if( split_diffusion ) { uK = KKn - ((xKL - xK) * lambdaK / dK); T = lambdaK / dK; } else //un-splitted diffusion { uK = KKn; T = 0; } real bcconds[3] = {0.0,1.0,0.0}; //pure neumann boundary condition if( tag_BC.isValid() && fKL.HaveData(tag_BC) ) //are there boundary conditions on face? { //retrive boundary conditions real_array bc = fKL.RealArray(tag_BC); bcconds[0] = bc[0]; bcconds[1] = bc[1]; bcconds[2] = bc[2]; } //account for boundary conditions R = A*T*bcconds[2]/(bcconds[0] + bcconds[1]*T + 1.0e-100); uK *=A* bcconds[0]/(bcconds[0] + bcconds[1]*T + 1.0e-100); T *=A* bcconds[0]/(bcconds[0] + bcconds[1]*T + 1.0e-100); //on BC if( uK.FrobeniusNorm() > 0 ) flag_DIFF = 1; else flag_DIFF = 0; } else std::cout << "No adjacent cell on non-boundary face" << std::endl; //record data for diffusion part if( tag_K.isValid() ) { fKL.Bulk(tag_DIFF_F) = flag_DIFF; fKL.Real(tag_DIFF_RT) = R; fKL.Real(tag_DIFF_CT) = T; real_array VT = fKL.RealArray(tag_DIFF_VT); VT[0] = uK(0,0); VT[1] = uK(0,1); VT[2] = uK(0,2); VT[3] = -uL(0,0); VT[4] = -uL(0,1); VT[5] = -uL(0,2); } //Advection part uK.Zero(); uL.Zero(); C = 0.0; cU = InvalidCell(); R = 0; if( U > 0.0 ) //flow out of back cell to front cell { cU = cK; C = U*A; //upstream corrector uK = (xK - xKL)*U*A; //downstream corrector if( cL.isValid() ) //internal face { r0 = r = xK - xL; if( tag_K.isValid() ) //heterogeneous media { KD = KL - KK; if( KD.FrobeniusNorm() > 0.0 ) { KD /= lambdaK + eps; iT = (rMatrix::Unit(3) + nKL.Transpose()*nKL*KD); iC = -nKL*dL*KD; r = (r*iT - iC)*(lambdaK/(lambdaK+eps)) + r*(eps/(lambdaK+eps)); //r = (r*iT - iC)*(lambdaK/(lambdaK+eps)) + (r + nKL*(KL-KK)/(U+eps))*(eps/(lambdaK+eps)); } } uL = (xKL - xL - r)*U*A; flag_CONV = 2; //internal } else if( fKL.GetMarker(boundary_face) ) flag_CONV = 1; //on outlet BC else std::cout << "No adjacent cell on non-boundary face" << std::endl; } else if( U < 0.0 ) //flow out of front cell to back cell { if( cL.isValid() ) //internal face { cU = cL; C = U*A; uL = (xKL - xL)*U*A; r = xL - xK; if( tag_K.isValid() ) //heterogeneous media { KD = KK - KL; if( KD.FrobeniusNorm() > 0.0 ) { KD /= lambdaL + eps; iT = (rMatrix::Unit(3) + nKL.Transpose()*nKL*KD); iC = nKL*dK*KD; r = (r*iT - iC)*(lambdaL/(lambdaL+eps)) + r*(eps/(lambdaL+eps)); //r = (r*iT - iC)*(lambdaL/(lambdaL+eps)) + (r + nKL*(KK-KL)/(U-eps))*(eps/(lambdaL+eps)); } } uK = (r + xK - xKL)*U*A; flag_CONV = 2; //internal } else if( fKL.GetMarker(boundary_face) )//boundary face { real bcconds[3] = {0.0,1.0,0.0}; //pure neumann boundary condition if( tag_BC.isValid() && fKL.HaveData(tag_BC) ) //are there boundary conditions on face? { //retrive boundary conditions real_array bc = fKL.RealArray(tag_BC); bcconds[0] = bc[0]; bcconds[1] = bc[1]; bcconds[2] = bc[2]; } if( bcconds[0] > 0 && fabs(bcconds[1]) < 1.0e-7 ) //Dirichlet BC { C = 0; R = U*A*bcconds[2]/bcconds[0]; flag_CONV = 0; //on inlet BC } else if( tag_K.isValid() ) //get expression out of diffusion { v = - A * (KKn - (xKL - xK) * lambdaK / dK); T = A * lambdaK / dK; cU = cK; C = U*A*T*bcconds[1]/(bcconds[0] + bcconds[1]*T + 1.0e-100); R = U*A *bcconds[2]/(bcconds[0] + bcconds[1]*T + 1.0e-100); uK =U*A*v*bcconds[1]/(bcconds[0] + bcconds[1]*T + 1.0e-100); flag_CONV = 1; //treat as outlet BC } else { if( bcconds[0] > 0 && bcconds[1] > 0 ) //Mixed BC { initialization_failure = true; std::cout << "Mixed BC type on inlet boundary face " << fKL->LocalID() << " is inconsistent for pure advection" << std::endl; } else if( bcconds[1] > 0 ) // Neumann BC { initialization_failure = true; std::cout << "Neumann BC type on inlet boundary face " << fKL->LocalID() << " is inconsistent for pure advection" << std::endl; } } } else std::cout << "No adjacent cell on non-boundary face" << std::endl; } else flag_CONV = 0; //otherwise velocity is zero, no flow through face //record data for advection part if( tag_U.isValid() ) { fKL.Bulk(tag_CONV_F) = flag_CONV; fKL.Real(tag_CONV_CU) = C; fKL.Reference(tag_CONV_EU) = cU.GetHandle(); fKL.Real(tag_CONV_RU) = R; real_array VU = fKL.RealArray(tag_CONV_VU); VU[0] = uK(0,0); VU[1] = uK(0,1); VU[2] = uK(0,2); VU[3] = uL(0,0); VU[4] = uL(0,1); VU[5] = uL(0,2); } // tag_U is defined } // for q } //openmp Tag failure = m->CreateTag("FAILURE",DATA_INTEGER,CELL,NONE,1); //Compute triplets integer negative = 0, total = 0; #if defined(USE_OMP) #pragma omp parallel reduction(+:negative,total) #endif { Cell cK; //current cell Face fKL; ElementArray<Face> faces; //faces of the cell std::vector<Stencil> compute; //approximations to be computed Tag tag_iC, tag_iT; //temporary data used for interpolation tensors { std::stringstream name; name << "INTERPOLATION_TENSOR_" << m->GetLocalProcessorRank(); tag_iT = m->CreateTag(name.str(),DATA_REAL,CELL,NONE); name.clear(); name << "INTERPOLATION_CORRECTION_" << m->GetLocalProcessorRank(); tag_iC = m->CreateTag(name.str(),DATA_REAL,CELL,NONE); } bulk flag; real_array v; #if defined(USE_OMP) #pragma omp for #endif for(integer q = 0; q < m->CellLastLocalID(); ++q ) if( m->isValidCell(q) ) { cK = m->CellByLocalID(q); compute.clear(); faces = cK.getFaces(); for(integer k = 0; k < faces.size(); ++k) { fKL = faces[k]; if( !BuildFlux(fKL) ) continue; //Diffusion stencils if( tag_K.isValid() ) { flag = fKL.Bulk(tag_DIFF_F); if( flag == 3 || flag == 1 ) { Stencil s; s.type = DIFFUSION; if( fKL.FaceOrientedOutside(cK) ) { s.coefs = fKL.RealArray(tag_DIFF_CB); s.elems = fKL.ReferenceArray(tag_DIFF_EB); s.rhs = &fKL.Real(tag_DIFF_RB); v = fKL.RealArray(tag_DIFF_VT); s.v[0] = v[0]; s.v[1] = v[1]; s.v[2] = v[2]; s.sign = 1; } else { s.coefs = fKL.RealArray(tag_DIFF_CF); s.elems = fKL.ReferenceArray(tag_DIFF_EF); s.rhs = &fKL.Real(tag_DIFF_RF); v = fKL.RealArray(tag_DIFF_VT); s.v[0] = v[3]; s.v[1] = v[4]; s.v[2] = v[5]; s.sign = -1; } if( s.v[0]*s.v[0] + s.v[1]*s.v[1] + s.v[2]*s.v[2] > 0.0 ) compute.push_back(s); } } //Advection stencils if( tag_U.isValid() ) { flag = fKL.Bulk(tag_CONV_F); if( flag >= 2 || flag == 1 ) { Stencil s; s.type = CONVECTION; if( fKL.FaceOrientedOutside(cK) ) { s.coefs = fKL.RealArray(tag_CONV_CB); s.elems = fKL.ReferenceArray(tag_CONV_EB); s.rhs = &fKL.Real(tag_CONV_RB); v = fKL.RealArray(tag_CONV_VU); s.v[0] = v[0]; s.v[1] = v[1]; s.v[2] = v[2]; s.sign = -1; } else { s.coefs = fKL.RealArray(tag_CONV_CF); s.elems = fKL.ReferenceArray(tag_CONV_EF); s.rhs = &fKL.Real(tag_CONV_RF); v = fKL.RealArray(tag_CONV_VU); s.v[0] = v[3]; s.v[1] = v[4]; s.v[2] = v[5]; s.sign = 1; } if( s.v[0]*s.v[0] + s.v[1]*s.v[1] + s.v[2]*s.v[2] > 0.0 ) compute.push_back(s); } } } //Gathered all the stencils at once if( !compute.empty() ) { if( !find_stencils(cK,compute,tag_BC,tag_K,tag_iT,tag_iC,boundary_face,bridge_layers,degenerate_diffusion_regularization,max_layers) ) { initialization_failure = true; cK.Integer(failure) = 1; std::cout << "Cannot find stencil for cell " << cK->LocalID() << std::endl; for(size_t it = 0; it < compute.size(); ++it) { if( !compute[it].computed ) { std::cout << (compute[it].type == DIFFUSION ? "Diffusion" : "Advection"); std::cout << " vec (" << compute[it].v[0] << "," << compute[it].v[1] << "," << compute[it].v[2] << ")"; std::cout << std::endl; } } std::cout << "Total: " << compute.size() << std::endl; find_stencils(cK,compute,tag_BC,tag_K,tag_iT,tag_iC,boundary_face,bridge_layers,degenerate_diffusion_regularization,max_layers); } else for(size_t it = 0; it < compute.size(); ++it) if( !compute[it].nonnegative ) negative++; total += (integer)compute.size(); } } //end of loop over cells //release interpolation tags m->DeleteTag(tag_iC); m->DeleteTag(tag_iT); } //openmp negative = m->Integrate(negative); total = m->Integrate(total); if( m->GetProcessorRank() == 0 ) std::cout << " negative flux approximation: " << negative << "/" << total << std::endl; { //Synchronize initialization_failure integer sync = (initialization_failure ? 1 : 0); sync = m->Integrate(sync); if( sync ) initialization_failure = true; } if( !initialization_failure ) m->DeleteTag(failure); }