diff --git a/src/ASM/ASMbase.h b/src/ASM/ASMbase.h index 95f98c99..07736746 100644 --- a/src/ASM/ASMbase.h +++ b/src/ASM/ASMbase.h @@ -367,9 +367,9 @@ public: virtual void initMADOF(const int*) {} //! \brief Generates element groups for multi-threading of interior integrals. - virtual void generateThreadGroups(const Integrand&, bool = false) {} + virtual void generateThreadGroups(const Integrand&, bool = false, bool = false) {} //! \brief Generates element groups for multi-threading of boundary integrals. - virtual void generateThreadGroups(char, bool = false) {} + virtual void generateThreadGroups(char, bool = false, bool = false) {} // Methods for integration of finite element quantities. diff --git a/src/ASM/ASMs2D.C b/src/ASM/ASMs2D.C index e165d8e9..2c3d0280 100644 --- a/src/ASM/ASMs2D.C +++ b/src/ASM/ASMs2D.C @@ -2591,16 +2591,18 @@ bool ASMs2D::evalSolution (Matrix& sField, const IntegrandBase& integrand, } -void ASMs2D::generateThreadGroups (const Integrand& integrand, bool silence) +void ASMs2D::generateThreadGroups (const Integrand& integrand, bool silence, + bool ignoreGlobalLM) { const int p1 = surf->order_u() - 1; const int p2 = surf->order_v() - 1; - generateThreadGroups(p1, p2, silence); + generateThreadGroups(p1, p2, silence, ignoreGlobalLM); } -void ASMs2D::generateThreadGroups (size_t strip1, size_t strip2, bool silence) +void ASMs2D::generateThreadGroups (size_t strip1, size_t strip2, + bool silence, bool ignoreGlobalLM) { const int n1 = surf->numCoefs_u(); const int n2 = surf->numCoefs_v(); @@ -2652,7 +2654,8 @@ void ASMs2D::generateThreadGroups (size_t strip1, size_t strip2, bool silence) // Verify that the nodes on this thread are not present on the others for (k = 0; k < j; k++) for (nit = nodes[j].begin(); nit != nodes[j].end(); ++nit) - if (nodes[k].find(*nit) != nodes[k].end()) + if ((this->getLMType(*nit+1) != 'G' || !ignoreGlobalLM) && + nodes[k].find(*nit) != nodes[k].end()) std::cout <<"\n ** Warning: Node "<< *nit <<" is present on both" <<" thread "<< k+1 <<" and thread "<< j+1; } diff --git a/src/ASM/ASMs2D.h b/src/ASM/ASMs2D.h index bf0fd75f..8ade565f 100644 --- a/src/ASM/ASMs2D.h +++ b/src/ASM/ASMs2D.h @@ -526,13 +526,17 @@ protected: //! \brief Generates element groups for multi-threading of interior integrals. //! \param[in] integrand Object with problem-specific data and methods //! \param[in] silence If \e true, suppress threading group outprint - virtual void generateThreadGroups(const Integrand& integrand, bool silence); + //! \param[in] ignoreGlobalLM If \e true ignore global multipliers in sanity check + virtual void generateThreadGroups(const Integrand& integrand, bool silence, + bool ignoreGlobalLM); //! \brief Generates element groups for multi-threading of interior integrals. //! \param[in] strip1 Strip width in first direction //! \param[in] strip2 Strip width in second direction //! \param[in] silence If \e true, suppress threading group outprint - void generateThreadGroups(size_t strip1, size_t strip2, bool silence); + //! \param[in] ignoreGlobalLM If \e true ignore global multipliers in sanity check + void generateThreadGroups(size_t strip1, size_t strip2, + bool silence, bool ignoreGlobalLM); public: //! \brief Auxilliary function for computation of basis function indices. diff --git a/src/ASM/ASMs2Dmx.C b/src/ASM/ASMs2Dmx.C index 2e69e0c4..74ff5e3a 100644 --- a/src/ASM/ASMs2Dmx.C +++ b/src/ASM/ASMs2Dmx.C @@ -999,7 +999,8 @@ bool ASMs2Dmx::evalSolution (Matrix& sField, const IntegrandBase& integrand, } -void ASMs2Dmx::generateThreadGroups (const Integrand& integrand, bool silence) +void ASMs2Dmx::generateThreadGroups (const Integrand& integrand, bool silence, + bool ignoreGlobalLM) { int p1 = 0, p2 = 0; for (const auto& it : m_basis) { @@ -1009,7 +1010,7 @@ void ASMs2Dmx::generateThreadGroups (const Integrand& integrand, bool silence) p2 = it->order_v(); } - ASMs2D::generateThreadGroups(p1-1, p2-1, silence); + ASMs2D::generateThreadGroups(p1-1, p2-1, silence, ignoreGlobalLM); } diff --git a/src/ASM/ASMs2Dmx.h b/src/ASM/ASMs2Dmx.h index bf28e251..0e2beece 100644 --- a/src/ASM/ASMs2Dmx.h +++ b/src/ASM/ASMs2Dmx.h @@ -198,7 +198,9 @@ public: //! \brief Generates element groups for multi-threading of interior integrals. //! \param[in] integrand Object with problem-specific data and methods //! \param[in] silence If \e true, suppress threading group outprint - virtual void generateThreadGroups(const Integrand& integrand, bool silence); + //! \param[in] ignoreGlobalLM If \e true, ignore global multipliers in sanity check + virtual void generateThreadGroups(const Integrand& integrand, bool silence, + bool ignoreGlobalLM); //! \brief Returns the number of nodal points in each parameter direction. //! \param[out] n1 Number of nodes in first (u) direction diff --git a/src/ASM/ASMs3D.C b/src/ASM/ASMs3D.C index 58a0e0e2..896fcbe4 100644 --- a/src/ASM/ASMs3D.C +++ b/src/ASM/ASMs3D.C @@ -2973,18 +2973,19 @@ bool ASMs3D::evalSolution (Matrix& sField, const IntegrandBase& integrand, } -void ASMs3D::generateThreadGroups (const Integrand& integrand, bool silence) +void ASMs3D::generateThreadGroups (const Integrand& integrand, bool silence, + bool ignoreGlobalLM) { const int p1 = svol->order(0) - 1; const int p2 = svol->order(1) - 1; const int p3 = svol->order(2) - 1; - ASMs3D::generateThreadGroups(p1, p2, p3, silence); + ASMs3D::generateThreadGroups(p1, p2, p3, silence, ignoreGlobalLM); } void ASMs3D::generateThreadGroups (size_t strip1, size_t strip2, size_t strip3, - bool silence) + bool silence, bool ignoreGlobalLM) { const int p1 = svol->order(0) - 1; const int p2 = svol->order(1) - 1; @@ -3041,7 +3042,8 @@ void ASMs3D::generateThreadGroups (size_t strip1, size_t strip2, size_t strip3, // Verify that the nodes on this thread are not present on the others for (k = 0; k < j; k++) for (nit = nodes[j].begin(); nit != nodes[j].end(); ++nit) - if (nodes[k].find(*nit) != nodes[k].end()) + if ((this->getLMType(*nit+1) != 'G' || !ignoreGlobalLM) && + nodes[k].find(*nit) != nodes[k].end()) std::cout <<"\n ** Warning: Node "<< *nit <<" is present on both" <<" thread "<< k+1 <<" and thread "<< j+1; } @@ -3050,7 +3052,7 @@ void ASMs3D::generateThreadGroups (size_t strip1, size_t strip2, size_t strip3, } -void ASMs3D::generateThreadGroups (char lIndex, bool silence) +void ASMs3D::generateThreadGroups (char lIndex, bool silence, bool) { if (threadGroupsFace.find(lIndex) != threadGroupsFace.end()) return; diff --git a/src/ASM/ASMs3D.h b/src/ASM/ASMs3D.h index 12d5d856..c8ebef98 100644 --- a/src/ASM/ASMs3D.h +++ b/src/ASM/ASMs3D.h @@ -594,17 +594,21 @@ protected: //! \brief Generates element groups for multi-threading of interior integrals. //! \param[in] integrand Object with problem-specific data and methods //! \param[in] silence If \e true, suppress threading group outprint - virtual void generateThreadGroups(const Integrand& integrand, bool silence); + //! \param[in] ignoreGlobalLM If \e true, ignore global multipliers in sanity check + virtual void generateThreadGroups(const Integrand& integrand, bool silence, + bool ignoreGlobalLM); //! \brief Generates element groups for multi-threading of interior integrals. //! \param[in] strip1 Strip width in first direction //! \param[in] strip2 Strip width in second direction //! \param[in] strip3 Strip width in third direction //! \param[in] silence If \e true, suppress threading group outprint - void generateThreadGroups(size_t strip1, size_t strip2, size_t strip3, bool silence); + //! \param[in] ignoreGlobalLM If \e true, ignore global multipliers in sanity check + void generateThreadGroups(size_t strip1, size_t strip2, size_t strip3, + bool silence, bool ignoreGlobalLM); //! \brief Generates element groups for multi-threading of boundary integrals. //! \param[in] lIndex Local index [1,6] of the boundary face //! \param[in] silence If \e true, suppress threading group outprint - virtual void generateThreadGroups(char lIndex, bool silence); + virtual void generateThreadGroups(char lIndex, bool silence, bool); public: //! \brief Auxilliary function for computation of basis function indices. diff --git a/src/ASM/ASMs3Dmx.C b/src/ASM/ASMs3Dmx.C index a2ce8de1..465c87c0 100644 --- a/src/ASM/ASMs3Dmx.C +++ b/src/ASM/ASMs3Dmx.C @@ -1008,7 +1008,8 @@ bool ASMs3Dmx::evalSolution (Matrix& sField, const IntegrandBase& integrand, } -void ASMs3Dmx::generateThreadGroups (const Integrand& integrand, bool silence) +void ASMs3Dmx::generateThreadGroups (const Integrand& integrand, bool silence, + bool ignoreGlobalLM) { int p[3] = { 0, 0, 0 }; for (const auto& it : m_basis) @@ -1016,16 +1017,16 @@ void ASMs3Dmx::generateThreadGroups (const Integrand& integrand, bool silence) if (it->order(d) > p[d]) p[d] = it->order(d); - ASMs3D::generateThreadGroups(p[0]-1,p[1]-1,p[2]-1,silence); + ASMs3D::generateThreadGroups(p[0]-1,p[1]-1,p[2]-1,silence,ignoreGlobalLM); } -void ASMs3Dmx::generateThreadGroups (char lIndex, bool silence) +void ASMs3Dmx::generateThreadGroups (char lIndex, bool silence, bool) { #ifdef USE_OPENMP omp_set_num_threads(1); #endif - ASMs3D::generateThreadGroups(lIndex,silence); + ASMs3D::generateThreadGroups(lIndex,silence,false); } diff --git a/src/ASM/ASMs3Dmx.h b/src/ASM/ASMs3Dmx.h index 6e4a29b8..89829b2f 100644 --- a/src/ASM/ASMs3Dmx.h +++ b/src/ASM/ASMs3Dmx.h @@ -191,11 +191,13 @@ public: //! \brief Generates element groups for multi-threading of interior integrals. //! \param[in] integrand Object with problem-specific data and methods //! \param[in] silence If \e true, suppress threading group outprint - virtual void generateThreadGroups(const Integrand& integrand, bool silence); + //! \param[in] ignoreGlobalLM If \e true ignore global multipliers in sanity check + virtual void generateThreadGroups(const Integrand& integrand, bool silence, + bool ignoreGlobalLM); //! \brief Generates element groups for multi-threading of boundary integrals. //! \param[in] lIndex Local index [1,6] of the boundary face //! \param[in] silence If \e true, suppress threading group outprint - virtual void generateThreadGroups(char lIndex, bool silence); + virtual void generateThreadGroups(char lIndex, bool silence, bool); //! \brief Returns the number of nodal points in each parameter direction. //! \param[out] n1 Number of nodes in first (u) direction diff --git a/src/SIM/SIMbase.C b/src/SIM/SIMbase.C index 50afa3f1..75e39404 100644 --- a/src/SIM/SIMbase.C +++ b/src/SIM/SIMbase.C @@ -57,6 +57,7 @@ SIMbase::SIMbase (IntegrandBase* itg) : g2l(&myGlb2Loc) myGen = nullptr; nGlPatches = 0; nIntGP = nBouGP = 0; + lagMTOK = false; MPCLess::compareSlaveDofOnly = true; // to avoid multiple slave definitions } @@ -1066,7 +1067,7 @@ bool SIMbase::preprocess (const IntVec& ignored, bool fixDup) bool silence = msgLevel < 1 || (msgLevel < 2 && myModel.size() > 1); for (mit = myModel.begin(); mit != myModel.end() && myProblem; ++mit) if (!(*mit)->empty()) - (*mit)->generateThreadGroups(*myProblem,silence); + (*mit)->generateThreadGroups(*myProblem,silence,lagMTOK); for (q = myProps.begin(); q != myProps.end(); ++q) if (q->pcode == Property::NEUMANN || @@ -1112,7 +1113,7 @@ void SIMbase::generateThreadGroups (const Property& p, bool silence) { ASMbase* pch = this->getPatch(p.patch); if (pch && abs(p.ldim)+1 == pch->getNoParamDim()) - pch->generateThreadGroups(p.lindx,silence); + pch->generateThreadGroups(p.lindx,silence,lagMTOK); } diff --git a/src/SIM/SIMbase.h b/src/SIM/SIMbase.h index dc40e17a..e2c68f0c 100644 --- a/src/SIM/SIMbase.h +++ b/src/SIM/SIMbase.h @@ -740,6 +740,7 @@ protected: // Model attributes bool isRefined; //!< Indicates if the model is adaptively refined + bool lagMTOK; //!< Indicates that global multipliers is okay with multithreading (app specific). unsigned char nsd; //!< Number of spatial dimensions PatchVec myModel; //!< The actual NURBS/spline model PropertyVec myProps; //!< Physical property mapping