bool ASMs3DLag::integrate (Integrand& integrand, int lIndex, GlobalIntegral& glInt, const TimeDomain& time) { if (!svol) return true; // silently ignore empty patches std::map<char,ThreadGroups>::const_iterator tit; if ((tit = threadGroupsFace.find(lIndex)) == threadGroupsFace.end()) { std::cerr <<" *** ASMs3DLag::integrate: No thread groups for face "<< lIndex << std::endl; return false; } const ThreadGroups& threadGrp = tit->second; // Get Gaussian quadrature points and weights int nGP = integrand.getBouIntegrationPoints(nGauss); const double* xg = GaussQuadrature::getCoord(nGP); const double* wg = GaussQuadrature::getWeight(nGP); if (!xg || !wg) return false; // Find the parametric direction of the face normal {-3,-2,-1, 1, 2, 3} const int faceDir = (lIndex+1)/(lIndex%2 ? -2 : 2); const int t0 = abs(faceDir); // unsigned normal direction of the face const int t1 = 1 + t0%3; // first tangent direction of the face const int t2 = 1 + t1%3; // second tangent direction of the face // Order of basis in the three parametric directions (order = degree + 1) const int p1 = svol->order(0); const int p2 = svol->order(1); const int p3 = svol->order(2); // Number of elements in each direction const int nel1 = (nx-1)/(p1-1); const int nel2 = (ny-1)/(p2-1); const int nel3 = (nz-1)/(p3-1); // Get parametric coordinates of the elements RealArray upar, vpar, wpar; if (t0 == 1) upar.resize(1,faceDir < 0 ? svol->startparam(0) : svol->endparam(0)); else if (t0 == 2) vpar.resize(1,faceDir < 0 ? svol->startparam(1) : svol->endparam(1)); else if (t0 == 3) wpar.resize(1,faceDir < 0 ? svol->startparam(2) : svol->endparam(2)); if (upar.empty()) this->getGridParameters(upar,0,1); if (vpar.empty()) this->getGridParameters(vpar,1,1); if (wpar.empty()) this->getGridParameters(wpar,2,1); // Integrate the extraordinary elements? size_t doXelms = 0; if (integrand.getIntegrandType() & Integrand::XO_ELEMENTS) if ((doXelms = nel1*nel2*nel3)*2 > MNPC.size()) { std::cerr <<" *** ASMs2DLag::integrate: Too few XO-elements " << MNPC.size() - doXelms << std::endl; return false; } std::map<char,size_t>::const_iterator iit = firstBp.find(lIndex); size_t firstp = iit == firstBp.end() ? 0 : iit->second; // === Assembly loop over all elements on the patch face ===================== bool ok = true; for (size_t g = 0; g < threadGrp.size() && ok; g++) { #pragma omp parallel for schedule(static) for (size_t t = 0; t < threadGrp[g].size(); t++) { FiniteElement fe(p1*p2*p3); fe.u = upar.front(); fe.v = vpar.front(); fe.w = wpar.front(); Matrix dNdu, Xnod, Jac; Vec4 X; Vec3 normal; double xi[3]; for (size_t l = 0; l < threadGrp[g][t].size() && ok; l++) { int iel = threadGrp[g][t][l]; int i1 = iel % nel1; int i2 = (iel / nel1) % nel2; int i3 = iel / (nel1*nel2); // Set up nodal point coordinates for current element if (!this->getElementCoordinates(Xnod,++iel)) { ok = false; break; } // Initialize element quantities fe.iel = abs(MLGE[doXelms+iel-1]); LocalIntegral* A = integrand.getLocalIntegral(fe.N.size(),fe.iel,true); if (!integrand.initElementBou(MNPC[doXelms+iel-1],*A)) { A->destruct(); ok = false; break; } // Define some loop control variables depending on which face we are on int nf1, j1, j2; switch (abs(faceDir)) { case 1: nf1 = nel2; j2 = i3; j1 = i2; break; case 2: nf1 = nel1; j2 = i3; j1 = i1; break; case 3: nf1 = nel1; j2 = i2; j1 = i1; break; default: nf1 = j1 = j2 = 0; } // --- Integration loop over all Gauss points in each direction -------- int k1, k2, k3; int jp = (j2*nf1 + j1)*nGP*nGP; fe.iGP = firstp + jp; // Global integration point counter for (int j = 0; j < nGP; j++) for (int i = 0; i < nGP; i++, fe.iGP++) { // Local element coordinates of current integration point xi[t0-1] = faceDir < 0 ? -1.0 : 1.0; xi[t1-1] = xg[i]; xi[t2-1] = xg[j]; fe.xi = xi[0]; fe.eta = xi[1]; fe.zeta = xi[2]; // Local element coordinates and parameter values // of current integration point switch (abs(faceDir)) { case 1: k2 = i; k3 = j; k1 = -1; break; case 2: k1 = i; k3 = j; k2 = -1; break; case 3: k1 = i; k2 = j; k3 = -1; break; default: k1 = k2 = k3 = -1; } if (upar.size() > 1) fe.u = 0.5*(upar[i1]*(1.0-xg[k1]) + upar[i1+1]*(1.0+xg[k1])); if (vpar.size() > 1) fe.v = 0.5*(vpar[i2]*(1.0-xg[k2]) + vpar[i2+1]*(1.0+xg[k2])); if (wpar.size() > 1) fe.w = 0.5*(wpar[i3]*(1.0-xg[k3]) + wpar[i3+1]*(1.0+xg[k3])); // Compute the basis functions and their derivatives, using // tensor product of one-dimensional Lagrange polynomials if (!Lagrange::computeBasis(fe.N,dNdu,p1,xi[0],p2,xi[1],p3,xi[2])) ok = false; // Compute basis function derivatives and the face normal fe.detJxW = utl::Jacobian(Jac,normal,fe.dNdX,Xnod,dNdu,t1,t2); if (fe.detJxW == 0.0) continue; // skip singular points if (faceDir < 0) normal *= -1.0; // Cartesian coordinates of current integration point X = Xnod * fe.N; X.t = time.t; // Evaluate the integrand and accumulate element contributions fe.detJxW *= wg[i]*wg[j]; if (!integrand.evalBou(*A,fe,time,X,normal)) ok = false; } // Finalize the element quantities if (ok && !integrand.finalizeElementBou(*A,fe,time)) ok = false; // Assembly of global system integral if (ok && !glInt.assemble(A->ref(),fe.iel)) ok = false; A->destruct(); } } } return ok; }
bool ASMs2DLag::integrate (Integrand& integrand, int lIndex, GlobalIntegral& glInt, const TimeDomain& time) { if (this->empty()) return true; // silently ignore empty patches // Get Gaussian quadrature points and weights int nG1 = this->getNoGaussPt(lIndex%10 < 3 ? p1 : p2, true); int nGP = integrand.getBouIntegrationPoints(nG1); const double* xg = GaussQuadrature::getCoord(nGP); const double* wg = GaussQuadrature::getWeight(nGP); if (!xg || !wg) return false; // Find the parametric direction of the edge normal {-2,-1, 1, 2} const int edgeDir = (lIndex%10+1)/((lIndex%2) ? -2 : 2); const int t1 = abs(edgeDir); // tangent direction normal to the patch edge const int t2 = 3-t1; // tangent direction along the patch edge // Number of elements in each direction const int nelx = (nx-1)/(p1-1); const int nely = (ny-1)/(p2-1); // Get parametric coordinates of the elements FiniteElement fe(p1*p2); RealArray upar, vpar; if (t1 == 1) { fe.u = edgeDir < 0 ? surf->startparam_u() : surf->endparam_u(); this->getGridParameters(vpar,1,1); } else if (t1 == 2) { this->getGridParameters(upar,0,1); fe.v = edgeDir < 0 ? surf->startparam_v() : surf->endparam_v(); } // Extract the Neumann order flag (1 or higher) for the integrand integrand.setNeumannOrder(1 + lIndex/10); // Integrate the extraordinary elements? size_t doXelms = 0; if (integrand.getIntegrandType() & Integrand::XO_ELEMENTS) if ((doXelms = nelx*nely)*2 > MNPC.size()) { std::cerr <<" *** ASMs2DLag::integrate: Too few XO-elements " << MNPC.size() - doXelms << std::endl; return false; } std::map<char,size_t>::const_iterator iit = firstBp.find(lIndex%10); size_t firstp = iit == firstBp.end() ? 0 : iit->second; Matrix dNdu, Xnod, Jac; Vec4 X; Vec3 normal; double xi[2]; // === Assembly loop over all elements on the patch edge ===================== int iel = 1; for (int i2 = 0; i2 < nely; i2++) for (int i1 = 0; i1 < nelx; i1++, iel++) { // Skip elements that are not on current boundary edge bool skipMe = false; switch (edgeDir) { case -1: if (i1 > 0) skipMe = true; break; case 1: if (i1 < nelx-1) skipMe = true; break; case -2: if (i2 > 0) skipMe = true; break; case 2: if (i2 < nely-1) skipMe = true; break; } if (skipMe) continue; // Set up nodal point coordinates for current element if (!this->getElementCoordinates(Xnod,iel)) return false; // Initialize element quantities fe.iel = abs(MLGE[doXelms+iel-1]); LocalIntegral* A = integrand.getLocalIntegral(fe.N.size(),fe.iel,true); bool ok = integrand.initElementBou(MNPC[doXelms+iel-1],*A); // --- Integration loop over all Gauss points along the edge ------------- int jp = (t1 == 1 ? i2 : i1)*nGP; fe.iGP = firstp + jp; // Global integration point counter for (int i = 0; i < nGP && ok; i++, fe.iGP++) { // Local element coordinates of current integration point xi[t1-1] = edgeDir < 0 ? -1.0 : 1.0; xi[t2-1] = xg[i]; fe.xi = xi[0]; fe.eta = xi[1]; // Parameter values of current integration point if (upar.size() > 1) fe.u = 0.5*(upar[i1]*(1.0-xg[i]) + upar[i1+1]*(1.0+xg[i])); if (vpar.size() > 1) fe.v = 0.5*(vpar[i2]*(1.0-xg[i]) + vpar[i2+1]*(1.0+xg[i])); // Compute the basis functions and their derivatives, using // tensor product of one-dimensional Lagrange polynomials if (!Lagrange::computeBasis(fe.N,dNdu,p1,xi[0],p2,xi[1])) ok = false; // Compute basis function derivatives and the edge normal fe.detJxW = utl::Jacobian(Jac,normal,fe.dNdX,Xnod,dNdu,t1,t2); if (fe.detJxW == 0.0) continue; // skip singular points if (edgeDir < 0) normal *= -1.0; // Cartesian coordinates of current integration point X = Xnod * fe.N; X.t = time.t; // Evaluate the integrand and accumulate element contributions fe.detJxW *= wg[i]; if (ok && !integrand.evalBou(*A,fe,time,X,normal)) ok = false; } // Finalize the element quantities if (ok && !integrand.finalizeElementBou(*A,fe,time)) ok = false; // Assembly of global system integral if (ok && !glInt.assemble(A->ref(),fe.iel)) ok = false; A->destruct(); if (!ok) return false; } return true; }