Skip to content

Commit

Permalink
adoption of new domain iterators
Browse files Browse the repository at this point in the history
  • Loading branch information
hverhelst committed Feb 14, 2025
1 parent fe86c6f commit a5040cb
Show file tree
Hide file tree
Showing 8 changed files with 68 additions and 51 deletions.
20 changes: 11 additions & 9 deletions src/gsC1SurfBasisEdge.h
Original file line number Diff line number Diff line change
Expand Up @@ -263,22 +263,24 @@ namespace gismo
const gsGeometry<T> & patch = m_mp.patch(0);

// Initialize domain element iterator
typename gsBasis<T>::domainIter domIt = m_geo.basis(0).makeDomainIterator(boundary::none);

#ifdef _OPENMP
for ( domIt->next(tid); domIt->good(); domIt->next(nt) )
#else
for (; domIt->good(); domIt->next() )
#endif
typename gsBasis<T>::domainIter domIt = m_geo.basis(0).domain()->beginBdr(boundary::none);
typename gsBasis<T>::domainIter domItEnd = m_geo.basis(0).domain()->endBdr(boundary::none);

# ifdef _OPENMP
domIt += tid;
for ( ; domIt<domItEnd; domIt+=nt )
# else
for (; domIt<domItEnd; ++domIt )
# endif
{
// Map the Quadrature rule to the element
quRule.mapTo( domIt->lowerCorner(), domIt->upperCorner(), quNodes, quWeights );
quRule.mapTo( domIt.lowerCorner(), domIt.upperCorner(), quNodes, quWeights );

// Perform required evaluations on the quadrature nodes
visitor_.evaluate(bf_index, typeBf, basis_g1, basis_geo, basis_plus, basis_minus, patch, quNodes, m_uv, m_gD, m_isBoundary);

// Assemble on element
visitor_.assemble(*domIt, quWeights);
visitor_.assemble(domIt, quWeights);

// Push to global matrix and right-hand side vector
#pragma omp critical(localToGlobal)
Expand Down
20 changes: 11 additions & 9 deletions src/gsC1SurfBasisVertex.h
Original file line number Diff line number Diff line change
Expand Up @@ -230,22 +230,24 @@ namespace gismo
const gsGeometry<T> & patch = m_mp.patch(0);

// Initialize domain element iterator
typename gsBasis<T>::domainIter domIt = m_geo.basis(0).makeDomainIterator(boundary::none);

#ifdef _OPENMP
for ( domIt->next(tid); domIt->good(); domIt->next(nt) )
#else
for (; domIt->good(); domIt->next() )
#endif
typename gsBasis<T>::domainIter domIt = m_geo.basis(0).domain()->beginBdr(boundary::none);
typename gsBasis<T>::domainIter domItEnd = m_geo.basis(0).domain()->endBdr(boundary::none);

# ifdef _OPENMP
domIt += tid;
for ( ; domIt<domItEnd; domIt+=nt )
# else
for (; domIt<domItEnd; ++domIt )
# endif
{
// Map the Quadrature rule to the element
quRule.mapTo( domIt->lowerCorner(), domIt->upperCorner(), quNodes, quWeights );
quRule.mapTo( domIt.lowerCorner(), domIt.upperCorner(), quNodes, quWeights );
#pragma omp critical(evaluate)
// Perform required evaluations on the quadrature nodes
visitor_.evaluate(basis_g1, basis_geo, m_basis_plus, m_basis_minus, patch, quNodes, m_gD, m_isBoundary, m_Phi);

// Assemble on element
visitor_.assemble(*domIt, quWeights);
visitor_.assemble(domIt, quWeights);

// Push to global matrix and right-hand side vector
#pragma omp critical(localToGlobal)
Expand Down
42 changes: 23 additions & 19 deletions src/gsC1SurfGluingData.h
Original file line number Diff line number Diff line change
Expand Up @@ -176,23 +176,25 @@ class gsC1SurfGluingData : public gsC1SurfGD<T>
// Initialize reference quadrature rule and visitor data
visitor_.initialize(basis,quRule);

// Initialize domain element iterator -- using unknown 0
typename gsBasis<T>::domainIter domIt = basis.makeDomainIterator(boundary::none);

#ifdef _OPENMP
for ( domIt->next(tid); domIt->good(); domIt->next(nt) )
#else
for (; domIt->good(); domIt->next() )
#endif
// Initialize domain element iterator
typename gsBasis<T>::domainIter domIt = basis.domain()->beginBdr(boundary::none);
typename gsBasis<T>::domainIter domItEnd = basis.domain()->endBdr(boundary::none);

# ifdef _OPENMP
domIt += tid;
for ( ; domIt<domItEnd; domIt+=nt )
# else
for (; domIt<domItEnd; ++domIt )
# endif
{
// Map the Quadrature rule to the element
quRule.mapTo( domIt->lowerCorner(), domIt->upperCorner(), quNodes, quWeights );
quRule.mapTo( domIt.lowerCorner(), domIt.upperCorner(), quNodes, quWeights );

// Perform required evaluations on the quadrature nodes
visitor_.evaluate(quNodes, this->m_mp);

// Assemble on element
visitor_.assemble(*domIt, quWeights);
visitor_.assemble(domIt, quWeights);

// Push to global matrix and right-hand side vector
#pragma omp critical(localToGlobal)
Expand Down Expand Up @@ -223,23 +225,25 @@ class gsC1SurfGluingData : public gsC1SurfGD<T>
// Initialize reference quadrature rule and visitor data
visitor_Beta.initialize(basis,quRule);

// Initialize domain element iterator -- using unknown 0
typename gsBasis<T>::domainIter domIt = basis.makeDomainIterator(boundary::none);
// Initialize domain element iterator
typename gsBasis<T>::domainIter domIt = basis.domain()->beginBdr(boundary::none);
typename gsBasis<T>::domainIter domItEnd = basis.domain()->endBdr(boundary::none);

#ifdef _OPENMP
for ( domIt->next(tid); domIt->good(); domIt->next(nt) )
#else
for (; domIt->good(); domIt->next() )
#endif
# ifdef _OPENMP
domIt += tid;
for ( ; domIt<domItEnd; domIt+=nt )
# else
for (; domIt<domItEnd; ++domIt )
# endif
{
// Map the Quadrature rule to the element
quRule.mapTo( domIt->lowerCorner(), domIt->upperCorner(), quNodes, quWeights );
quRule.mapTo( domIt.lowerCorner(), domIt.upperCorner(), quNodes, quWeights );

// Perform required evaluations on the quadrature nodes
visitor_Beta.evaluateBeta(quNodes, this->m_mp, sol);

// Assemble on element
visitor_Beta.assembleBeta(*domIt, quWeights);
visitor_Beta.assembleBeta(domIt, quWeights);

// Push to global matrix and right-hand side vector
#pragma omp critical(localToGlobal)
Expand Down
22 changes: 12 additions & 10 deletions src/gsC1SurfGluingDataAssembler.h
Original file line number Diff line number Diff line change
Expand Up @@ -127,7 +127,7 @@ class gsC1SurfGluingDataAssembler : public gsAssembler<T>
m_system_beta_S.matrix().makeCompressed();

}

template <class T, class bhVisitor>
inline void gsC1SurfGluingDataAssembler<T, bhVisitor>::apply(bhVisitor & visitor,
index_t patchIndex,
Expand Down Expand Up @@ -155,23 +155,25 @@ class gsC1SurfGluingDataAssembler : public gsAssembler<T>

//const gsGeometry<T> & patch = m_geo.patch(patchIndex); // 0 = patchindex

// Initialize domain element iterator -- using unknown 0
typename gsBasis<T>::domainIter domIt = basis.makeDomainIterator(boundary::none);
// Initialize domain element iterator
typename gsBasis<T>::domainIter domIt = basis.domain()->beginBdr(boundary::none);
typename gsBasis<T>::domainIter domItEnd = basis.domain()->endBdr(boundary::none);

#ifdef _OPENMP
for ( domIt->next(tid); domIt->good(); domIt->next(nt) )
#else
for (; domIt->good(); domIt->next() )
#endif
# ifdef _OPENMP
domIt += tid;
for ( ; domIt<domItEnd; domIt+=nt )
# else
for (; domIt<domItEnd; ++domIt )
# endif
{
// Map the Quadrature rule to the element
quRule.mapTo( domIt->lowerCorner(), domIt->upperCorner(), quNodes, quWeights );
quRule.mapTo( domIt.lowerCorner(), domIt.upperCorner(), quNodes, quWeights );

// Perform required evaluations on the quadrature nodes
visitor_.evaluate(basis, quNodes, m_uv, m_mp, m_gamma, m_isBoundary);

// Assemble on element
visitor_.assemble(*domIt, quWeights);
visitor_.assemble(domIt, quWeights);

// Push to global matrix and right-hand side vector
#pragma omp critical(localToGlobal)
Expand Down
4 changes: 2 additions & 2 deletions src/gsC1SurfGluingDataVisitor.h
Original file line number Diff line number Diff line change
Expand Up @@ -263,7 +263,7 @@ class gsC1SurfGluingDataVisitor
localRhsBeta.setZero(numActiveBeta, 1 ); //multiple right-hand sides
}

inline void assemble(gsDomainIterator<T> & element,
inline void assemble(gsDomainIteratorWrapper<T> & element,
const gsVector<T> & quWeights)
{
GISMO_UNUSED(element);
Expand All @@ -285,7 +285,7 @@ class gsC1SurfGluingDataVisitor
}
}

inline void assembleBeta(gsDomainIterator<T> & element,
inline void assembleBeta(gsDomainIteratorWrapper<T> & element,
const gsVector<T> & quWeights)
{
GISMO_UNUSED(element);
Expand Down
2 changes: 1 addition & 1 deletion src/gsC1SurfVisitorBasisEdge.h
Original file line number Diff line number Diff line change
Expand Up @@ -158,7 +158,7 @@ namespace gismo
} // Patch 1
} // evaluate1

inline void assemble(gsDomainIterator<T> & element,
inline void assemble(gsDomainIteratorWrapper<T> & element,
const gsVector<T> & quWeights)
{
GISMO_UNUSED(element);
Expand Down
2 changes: 1 addition & 1 deletion src/gsC1SurfVisitorBasisVertex.h
Original file line number Diff line number Diff line change
Expand Up @@ -419,7 +419,7 @@ namespace gismo

} // evaluate

inline void assemble(gsDomainIterator<T> & element,
inline void assemble(gsDomainIteratorWrapper<T> & element,
const gsVector<T> & quWeights)
{
GISMO_UNUSED(element);
Expand Down
7 changes: 7 additions & 0 deletions src/gsContainerBasis.h
Original file line number Diff line number Diff line change
Expand Up @@ -130,12 +130,19 @@ namespace gismo

index_t nPieces() const {return basisContainer.size();}

virtual memory::shared_ptr<gsDomain<T> > domain() const
{
return basisContainer[0].domain();
}

GISMO_DEPRECATED
typename gsBasis<T>::domainIter makeDomainIterator(const boxSide & side) const override
{
// Using the inner basis for iterating
return basisContainer[0].makeDomainIterator(side);
}

GISMO_DEPRECATED
typename gsBasis<T>::domainIter makeDomainIterator() const override
{
// Using the inner basis for iterating
Expand Down

0 comments on commit a5040cb

Please sign in to comment.