From 53f5d1f830653fe4365ca0f9548fbd7fad3bd168 Mon Sep 17 00:00:00 2001 From: Arne Morten Kvarving Date: Fri, 25 Aug 2023 10:48:53 +0200 Subject: [PATCH] changed: no reason to override write in ASMsxDmx we can just use the ASMsxD implementation --- src/ASM/ASMs2D.C | 8 +++++--- src/ASM/ASMs2Dmx.C | 15 --------------- src/ASM/ASMs2Dmx.h | 2 -- src/ASM/ASMs3D.C | 8 +++++--- src/ASM/ASMs3Dmx.C | 15 --------------- src/ASM/ASMs3Dmx.h | 2 -- 6 files changed, 10 insertions(+), 40 deletions(-) diff --git a/src/ASM/ASMs2D.C b/src/ASM/ASMs2D.C index bc7337a2..dfe0fd93 100644 --- a/src/ASM/ASMs2D.C +++ b/src/ASM/ASMs2D.C @@ -181,12 +181,14 @@ bool ASMs2D::read (std::istream& is) } -bool ASMs2D::write (std::ostream& os, int) const +bool ASMs2D::write (std::ostream& os, int basis) const { if (!surf) return false; + if (basis > static_cast(this->getNoBasis())) return false; + const Go::SplineSurface* spline = this->getBasis(basis); + if (!spline) return false; - os <<"200 1 0 0\n"; - os << *surf; + os <<"200 1 0 0\n" << *spline; return os.good(); } diff --git a/src/ASM/ASMs2Dmx.C b/src/ASM/ASMs2Dmx.C index a6adc81b..16e2b938 100644 --- a/src/ASM/ASMs2Dmx.C +++ b/src/ASM/ASMs2Dmx.C @@ -112,21 +112,6 @@ bool ASMs2Dmx::readBasis (std::istream& is, size_t basis) } -bool ASMs2Dmx::write (std::ostream& os, int basis) const -{ - if (basis == -1) - os <<"200 1 0 0\n" << *projB; - else if (basis < 1 || basis > (int)m_basis.size()) - os <<"200 1 0 0\n" << *surf; - else if (m_basis[basis-1]) - os <<"200 1 0 0\n" << *m_basis[basis-1]; - else - return false; - - return os.good(); -} - - void ASMs2Dmx::clear (bool retainGeometry) { // these are managed by shared ptrs, make sure base class do not delete them. diff --git a/src/ASM/ASMs2Dmx.h b/src/ASM/ASMs2Dmx.h index 2d09a978..1395c4bb 100644 --- a/src/ASM/ASMs2Dmx.h +++ b/src/ASM/ASMs2Dmx.h @@ -55,8 +55,6 @@ public: //! \brief Reads a basis from the given input stream. virtual bool readBasis(std::istream& is, size_t basis); - //! \brief Writes the geometry/basis of the patch to given stream. - virtual bool write(std::ostream& os, int basis) const; //! \brief Generates the finite element topology data for the patch. //! \details The data generated are the element-to-node connectivity array, diff --git a/src/ASM/ASMs3D.C b/src/ASM/ASMs3D.C index 5884b43b..6e4c04ab 100644 --- a/src/ASM/ASMs3D.C +++ b/src/ASM/ASMs3D.C @@ -152,12 +152,14 @@ bool ASMs3D::read (std::istream& is) } -bool ASMs3D::write (std::ostream& os, int) const +bool ASMs3D::write (std::ostream& os, int basis) const { if (!svol) return false; + if (basis > static_cast(this->getNoBasis())) return false; + const Go::SplineVolume* spline = this->getBasis(basis); + if (!spline) return false; - os <<"700 1 0 0\n"; - os << *svol; + os <<"700 1 0 0\n" << *spline; return os.good(); } diff --git a/src/ASM/ASMs3Dmx.C b/src/ASM/ASMs3Dmx.C index 5c3957b8..cbdba8ef 100644 --- a/src/ASM/ASMs3Dmx.C +++ b/src/ASM/ASMs3Dmx.C @@ -110,21 +110,6 @@ bool ASMs3Dmx::readBasis (std::istream& is, size_t basis) } -bool ASMs3Dmx::write (std::ostream& os, int basis) const -{ - if (basis == -1) - os <<"700 1 0 0\n" << *projB; - else if (basis < 1 || basis > (int)m_basis.size()) - os <<"700 1 0 0\n" << *svol; - else if (m_basis[basis-1]) - os <<"700 1 0 0\n" << *m_basis[basis-1]; - else - return false; - - return os.good(); -} - - void ASMs3Dmx::clear (bool retainGeometry) { // these are managed by shared ptrs, make sure base class do not delete them. diff --git a/src/ASM/ASMs3Dmx.h b/src/ASM/ASMs3Dmx.h index 3dc8bf9c..0ceac55c 100644 --- a/src/ASM/ASMs3Dmx.h +++ b/src/ASM/ASMs3Dmx.h @@ -74,8 +74,6 @@ public: //! \brief Reads a basis from the given input stream. virtual bool readBasis(std::istream& is, size_t basis); - //! \brief Writes the geometry/basis of the patch to given stream. - virtual bool write(std::ostream& os, int basis) const; //! \brief Returns the number of bases. virtual size_t getNoBasis() const { return m_basis.size(); }