From: geay Date: Tue, 20 May 2014 15:35:27 +0000 (+0200) Subject: On the road of the refinement for AMR mesh driven by a vector of bool. X-Git-Tag: V7_5_0a1~37^2~41 X-Git-Url: http://git.salome-platform.org/gitweb/?a=commitdiff_plain;h=36d31be36964def77645e16e439dcd008eb62339;p=modules%2Fmed.git On the road of the refinement for AMR mesh driven by a vector of bool. --- diff --git a/src/INTERP_KERNEL/BoxSplittingOptions.cxx b/src/INTERP_KERNEL/BoxSplittingOptions.cxx new file mode 100644 index 000000000..d932a482d --- /dev/null +++ b/src/INTERP_KERNEL/BoxSplittingOptions.cxx @@ -0,0 +1,45 @@ +// Copyright (C) 2007-2014 CEA/DEN, EDF R&D +// +// This library is free software; you can redistribute it and/or +// modify it under the terms of the GNU Lesser General Public +// License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// +// This library is distributed in the hope that it will be useful, +// but WITHOUT ANY WARRANTY; without even the implied warranty of +// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU +// Lesser General Public License for more details. +// +// You should have received a copy of the GNU Lesser General Public +// License along with this library; if not, write to the Free Software +// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA +// +// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com +// +// Author : Anthony Geay + +#include "BoxSplittingOptions.hxx" + +#include + +const double INTERP_KERNEL::BoxSplittingOptions::DFT_EFFECIENCY=0.5; + +const double INTERP_KERNEL::BoxSplittingOptions::DFT_EFFECIENCY_SND=0.5; + +void INTERP_KERNEL::BoxSplittingOptions::init() +{ + _effeciency=DFT_EFFECIENCY; + _effeciency_snd=DFT_EFFECIENCY_SND; + _min_cell_direction=DFT_MIN_CELL_DIRECTION; + _max_cells=DFT_MAX_CELLS; +} + +std::string INTERP_KERNEL::BoxSplittingOptions::printOptions() const +{ + std::ostringstream oss; + oss << "Efficiency threshold: " << 100*_effeciency << "%"<< std::endl; + oss << "Efficiency Snd threshold: " << 100*_effeciency_snd << "%"<< std::endl; + oss << "Min. box side length: " << _min_cell_direction << std::endl; + oss << "Max. box cells : " << _max_cells << std::endl; + return oss.str(); +} diff --git a/src/INTERP_KERNEL/BoxSplittingOptions.hxx b/src/INTERP_KERNEL/BoxSplittingOptions.hxx new file mode 100644 index 000000000..f8e47caed --- /dev/null +++ b/src/INTERP_KERNEL/BoxSplittingOptions.hxx @@ -0,0 +1,65 @@ +// Copyright (C) 2007-2014 CEA/DEN, EDF R&D +// +// This library is free software; you can redistribute it and/or +// modify it under the terms of the GNU Lesser General Public +// License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// +// This library is distributed in the hope that it will be useful, +// but WITHOUT ANY WARRANTY; without even the implied warranty of +// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU +// Lesser General Public License for more details. +// +// You should have received a copy of the GNU Lesser General Public +// License along with this library; if not, write to the Free Software +// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA +// +// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com +// +// Author : Anthony Geay + +#ifndef __BOXSPLITTINGOPTIONS_HXX__ +#define __BOXSPLITTINGOPTIONS_HXX__ + +#include "INTERPKERNELDefines.hxx" +#include "NormalizedUnstructuredMesh.hxx" + +#include + +namespace INTERP_KERNEL +{ + /*! + * \class BoxSplittingOptions + * Class defining the options for box splitting used for AMR algorithm like creation of patches following a criterion. + */ + class INTERPKERNEL_EXPORT BoxSplittingOptions + { + private: + double _effeciency; + double _effeciency_snd; + int _min_cell_direction; + int _max_cells; + public: + BoxSplittingOptions() { init(); } + void init(); + double getEffeciency() const { return _effeciency; } + void setEffeciency(double effeciency) { _effeciency=effeciency; } + double getEffeciencySnd() const { return _effeciency_snd; } + void setEffeciencySnd(double effeciencySnd) { _effeciency_snd=effeciencySnd; } + int getMinCellDirection() const { return _min_cell_direction; } + void setMinCellDirection(int minCellDirection) { _min_cell_direction=minCellDirection; } + int getMaxCells() const { return _max_cells; } + void setMaxCells(int maxCells) { _max_cells=maxCells; } + void copyOptions(const BoxSplittingOptions & other) { *this=other; } + std::string printOptions() const; + private: + static const int DFT_MIN_CELL_DIRECTION=3; + static const int DFT_MAX_CELLS=1000; + static const double DFT_EFFECIENCY; + static const double DFT_EFFECIENCY_SND; + public: + static const char PRECISION_STR[]; + }; +} + +#endif diff --git a/src/INTERP_KERNEL/CMakeLists.txt b/src/INTERP_KERNEL/CMakeLists.txt index fd36e53f3..5a87870d9 100644 --- a/src/INTERP_KERNEL/CMakeLists.txt +++ b/src/INTERP_KERNEL/CMakeLists.txt @@ -28,6 +28,7 @@ SET(interpkernel_SOURCES CellModel.cxx UnitTetraIntersectionBary.cxx InterpolationOptions.cxx + BoxSplittingOptions.cxx DirectedBoundingBox.cxx Interpolation2DCurve.cxx Interpolation3DSurf.cxx diff --git a/src/MEDCoupling/MEDCouplingCartesianAMRMesh.cxx b/src/MEDCoupling/MEDCouplingCartesianAMRMesh.cxx index bf06a977e..929adf7d6 100644 --- a/src/MEDCoupling/MEDCouplingCartesianAMRMesh.cxx +++ b/src/MEDCoupling/MEDCouplingCartesianAMRMesh.cxx @@ -33,16 +33,20 @@ using namespace ParaMEDMEM; * \param [in] mesh not null pointer of refined mesh replacing the cell range of \a father defined by the bottom left and top right just after. * \param [in] bottomLeftTopRight a vector equal to the space dimension of \a mesh that specifies for each dimension, the included cell start of the range for the first element of the pair, * a the end cell (\b excluded) of the range for the second element of the pair. + * \param [in] factors The refinement per axis relative to the father of \a this. */ -MEDCouplingCartesianAMRPatch::MEDCouplingCartesianAMRPatch(MEDCouplingCartesianAMRMesh *mesh, const std::vector< std::pair >& bottomLeftTopRight) +MEDCouplingCartesianAMRPatch::MEDCouplingCartesianAMRPatch(MEDCouplingCartesianAMRMesh *mesh, const std::vector< std::pair >& bottomLeftTopRight, const std::vector& factors) { if(!mesh) throw INTERP_KERNEL::Exception("EDCouplingCartesianAMRPatch constructor : input mesh is NULL !"); _mesh=mesh; _mesh->incrRef(); - int dim((int)bottomLeftTopRight.size()); - if(dim!=_mesh->getSpaceDimension()) + int dim((int)bottomLeftTopRight.size()),dimExp(_mesh->getSpaceDimension()); + if(dim!=dimExp) throw INTERP_KERNEL::Exception("MEDCouplingCartesianAMRPatch constructor : space dimension of father and input bottomLeft/topRight size mismatches !"); _bl_tr=bottomLeftTopRight; + if((int)factors.size()!=dimExp) + throw INTERP_KERNEL::Exception("MEDCouplingCartesianAMRPatch constructor : space dimension of father and input factors per axis size mismatches !"); + _factors=factors; } int MEDCouplingCartesianAMRPatch::getNumberOfCellsRecursiveWithOverlap() const @@ -60,9 +64,9 @@ int MEDCouplingCartesianAMRPatch::getMaxNumberOfLevelsRelativeToThis() const return _mesh->getMaxNumberOfLevelsRelativeToThis(); } -void MEDCouplingCartesianAMRPatch::addPatch(const std::vector< std::pair >& bottomLeftTopRight, int factor) +void MEDCouplingCartesianAMRPatch::addPatch(const std::vector< std::pair >& bottomLeftTopRight, const std::vector& factors) { - return _mesh->addPatch(bottomLeftTopRight,factor); + return _mesh->addPatch(bottomLeftTopRight,factors); } int MEDCouplingCartesianAMRPatch::getNumberOfOverlapedCellsForFather() const @@ -154,17 +158,142 @@ void MEDCouplingCartesianAMRMesh::detachFromFather() /*! * \param [in] bottomLeftTopRight a vector equal to the space dimension of \a mesh that specifies for each dimension, the included cell start of the range for the first element of the pair, * a the end cell (\b excluded) of the range for the second element of the pair. - * \param [in] factor The != 0 factor of refinement. + * \param [in] factors The != 0 factor of refinement per axis. */ -void MEDCouplingCartesianAMRMesh::addPatch(const std::vector< std::pair >& bottomLeftTopRight, int factor) +void MEDCouplingCartesianAMRMesh::addPatch(const std::vector< std::pair >& bottomLeftTopRight, const std::vector& factors) { MEDCouplingAutoRefCountObjectPtr mesh(static_cast(_mesh->buildStructuredSubPart(bottomLeftTopRight))); - mesh->refineWithFactor(factor); + mesh->refineWithFactor(factors); MEDCouplingAutoRefCountObjectPtr zeMesh(new MEDCouplingCartesianAMRMesh(this,mesh)); - MEDCouplingAutoRefCountObjectPtr elt(new MEDCouplingCartesianAMRPatch(zeMesh,bottomLeftTopRight)); + MEDCouplingAutoRefCountObjectPtr elt(new MEDCouplingCartesianAMRPatch(zeMesh,bottomLeftTopRight,factors)); _patches.push_back(elt); } +/// @cond INTERNAL + +class InternalPatch : public RefCountObjectOnly +{ +public: + InternalPatch():_nb_of_true(0) { } + double getEfficiency() const { return (double)_nb_of_true/(double)_crit.size(); } + int getNumberOfCells() const { return (int)_crit.size(); } + void setNumberOfTrue(int nboft) { _nb_of_true=nboft; } + std::vector& getCriterion() { return _crit; } + std::vector< std::pair >& getPart() { return _part; } + const std::vector< std::pair >& getConstPart() const { return _part; } + bool presenceOfTrue() const { return _nb_of_true>0; } +protected: + ~InternalPatch() { } +private: + int _nb_of_true; + std::vector _crit; + std::vector< std::pair > _part; +}; + +#if 0 +void dissectBigPatch (const Mesh& mesh, const Field& fieldFlag, const unsigned int minCellDirection, + const unsigned int big_dims, const int dissect_direction, int cut[3] ) const +{ + int cut_found = 0 ; + int cut_place = -1 ; + float * ratio = NULL ; + float * ratio_inner = NULL ; + + ratio = new float [big_dims-1]; + for(unsigned int id=0; id& factors) +{ + if(!criterion || !criterion->isAllocated()) + throw INTERP_KERNEL::Exception("MEDCouplingCartesianAMRMesh::createPatchesFromCriterion : the criterion DataArrayByte instance must be allocated and not NULL !"); + int nbCells(getNumberOfCellsAtCurrentLevel()); + if(nbCells!=criterion->getNumberOfTuples()) + throw INTERP_KERNEL::Exception("MEDCouplingCartesianAMRMesh::createPatchesFromCriterion : the number of tuples of criterion array must be equal to the number of cells at the current level !"); + _patches.clear(); + // + std::vector cgs(_mesh->getCellGridStructure()); + std::vector crit(criterion->toVectorOfBool());//check that criterion has one component. + std::vector< MEDCouplingAutoRefCountObjectPtr > listOfPatches,listOfPatchesOK; + // + MEDCouplingAutoRefCountObjectPtr p(new InternalPatch); + p->setNumberOfTrue(MEDCouplingStructuredMesh::FindMinimalPartOf(cgs,crit,p->getCriterion(),p->getPart())); + if(p->presenceOfTrue()) + listOfPatches.push_back(p); + while(!listOfPatches.empty()) + { + std::vector< MEDCouplingAutoRefCountObjectPtr > listOfPatchesTmp; + for(std::vector< MEDCouplingAutoRefCountObjectPtr >::iterator it=listOfPatches.begin();it!=listOfPatches.end();it++) + { + // + if((*it)->getEfficiency()>=bso.getEffeciency()) + { + if((*it)->getNumberOfCells()>=bso.getMaxCells()) + { + + } + } + } + listOfPatches=listOfPatchesTmp; + } + for(std::vector< MEDCouplingAutoRefCountObjectPtr >::const_iterator it=listOfPatchesOK.begin();it!=listOfPatchesOK.end();it++) + { + addPatch((*it)->getConstPart(),factors); + } +} + void MEDCouplingCartesianAMRMesh::removePatch(int patchId) { checkPatchId(patchId); diff --git a/src/MEDCoupling/MEDCouplingCartesianAMRMesh.hxx b/src/MEDCoupling/MEDCouplingCartesianAMRMesh.hxx index 3815681b6..5f6309f48 100644 --- a/src/MEDCoupling/MEDCouplingCartesianAMRMesh.hxx +++ b/src/MEDCoupling/MEDCouplingCartesianAMRMesh.hxx @@ -26,24 +26,26 @@ #include "MEDCouplingRefCountObject.hxx" #include "MEDCouplingAutoRefCountObjectPtr.hxx" +#include "BoxSplittingOptions.hxx" #include "InterpKernelException.hxx" namespace ParaMEDMEM { class MEDCouplingIMesh; class MEDCouplingUMesh; + class DataArrayByte; class MEDCouplingCartesianAMRMesh; /// @cond INTERNAL class MEDCouplingCartesianAMRPatch : public RefCountObject { public: - MEDCouplingCartesianAMRPatch(MEDCouplingCartesianAMRMesh *mesh, const std::vector< std::pair >& bottomLeftTopRight); + MEDCouplingCartesianAMRPatch(MEDCouplingCartesianAMRMesh *mesh, const std::vector< std::pair >& bottomLeftTopRight, const std::vector& factors); // direct forward to _mesh MEDCOUPLING_EXPORT int getNumberOfCellsRecursiveWithOverlap() const; MEDCOUPLING_EXPORT int getNumberOfCellsRecursiveWithoutOverlap() const; MEDCOUPLING_EXPORT int getMaxNumberOfLevelsRelativeToThis() const; - MEDCOUPLING_EXPORT void addPatch(const std::vector< std::pair >& bottomLeftTopRight, int factor); + MEDCOUPLING_EXPORT void addPatch(const std::vector< std::pair >& bottomLeftTopRight, const std::vector& factors); // end of direct forward to _mesh MEDCOUPLING_EXPORT int getNumberOfOverlapedCellsForFather() const; // basic set/get @@ -55,6 +57,7 @@ namespace ParaMEDMEM private: //! bottom left/top right cell range relative to \a _father std::vector< std::pair > _bl_tr; + std::vector _factors; MEDCouplingAutoRefCountObjectPtr _mesh; }; /// @endcond @@ -78,7 +81,8 @@ namespace ParaMEDMEM MEDCOUPLING_EXPORT const MEDCouplingCartesianAMRMesh *getFather() const; MEDCOUPLING_EXPORT const MEDCouplingCartesianAMRMesh *getGodFather() const; MEDCOUPLING_EXPORT void detachFromFather(); - MEDCOUPLING_EXPORT void addPatch(const std::vector< std::pair >& bottomLeftTopRight, int factor); + MEDCOUPLING_EXPORT void addPatch(const std::vector< std::pair >& bottomLeftTopRight, const std::vector& factors); + MEDCOUPLING_EXPORT void createPatchesFromCriterion(const INTERP_KERNEL::BoxSplittingOptions& bso, const DataArrayByte *criterion, const std::vector& factors); MEDCOUPLING_EXPORT void removePatch(int patchId); MEDCOUPLING_EXPORT int getNumberOfPatches() const; MEDCOUPLING_EXPORT const MEDCouplingCartesianAMRPatch *getPatch(int patchId) const; diff --git a/src/MEDCoupling/MEDCouplingIMesh.cxx b/src/MEDCoupling/MEDCouplingIMesh.cxx index dad9a4f13..142fac4b8 100644 --- a/src/MEDCoupling/MEDCouplingIMesh.cxx +++ b/src/MEDCoupling/MEDCouplingIMesh.cxx @@ -182,17 +182,27 @@ MEDCouplingCMesh *MEDCouplingIMesh::convertToCartesian() const * The origin of \a this will be not touched only spacing and node structure will be changed. * This method can be useful for AMR users. */ -void MEDCouplingIMesh::refineWithFactor(int factor) +void MEDCouplingIMesh::refineWithFactor(const std::vector& factors) { - if(factor==0) - throw INTERP_KERNEL::Exception("MEDCouplingIMesh::refineWithFactor : refinement factor must be != 0 !"); + if((int)factors.size()!=_space_dim) + throw INTERP_KERNEL::Exception("MEDCouplingIMesh::refineWithFactor : refinement factors must have size equal to spaceDim !"); checkCoherency(); - int factAbs(std::abs(factor)); - double fact2(1./(double)factor); - std::transform(_structure,_structure+_space_dim,_structure,std::bind2nd(std::plus(),-1)); - std::transform(_structure,_structure+_space_dim,_structure,std::bind2nd(std::multiplies(),factAbs)); - std::transform(_structure,_structure+_space_dim,_structure,std::bind2nd(std::plus(),1)); - std::transform(_dxyz,_dxyz+_space_dim,_dxyz,std::bind2nd(std::multiplies(),fact2)); + std::vector structure(_structure,_structure+3); + std::vector dxyz(_dxyz,_dxyz+3); + for(int i=0;i<_space_dim;i++) + { + if(factors[i]<=0) + { + std::ostringstream oss; oss << "MEDCouplingIMesh::refineWithFactor : factor for axis #" << i << " (" << factors[i] << ")is invalid ! Must be > 0 !"; + throw INTERP_KERNEL::Exception(oss.str().c_str()); + } + int factAbs(std::abs(factors[i])); + double fact2(1./(double)factors[i]); + structure[i]=(_structure[i]-1)*factAbs+1; + dxyz[i]=fact2*_dxyz[i]; + } + std::copy(structure.begin(),structure.end(),_structure); + std::copy(dxyz.begin(),dxyz.end(),_dxyz); declareAsNew(); } diff --git a/src/MEDCoupling/MEDCouplingIMesh.hxx b/src/MEDCoupling/MEDCouplingIMesh.hxx index 3a487a80e..896a7835d 100644 --- a/src/MEDCoupling/MEDCouplingIMesh.hxx +++ b/src/MEDCoupling/MEDCouplingIMesh.hxx @@ -46,7 +46,7 @@ namespace ParaMEDMEM MEDCOUPLING_EXPORT std::string getAxisUnit() const; MEDCOUPLING_EXPORT double getMeasureOfAnyCell() const; MEDCOUPLING_EXPORT MEDCouplingCMesh *convertToCartesian() const; - MEDCOUPLING_EXPORT void refineWithFactor(int factor); + MEDCOUPLING_EXPORT void refineWithFactor(const std::vector& factors); MEDCOUPLING_EXPORT static void CondenseFineToCoarse(DataArrayDouble *coarseDA, const std::vector& coarseSt, const DataArrayDouble *fineDA, const std::vector< std::pair >& fineLocInCoarse, const std::vector& facts); MEDCOUPLING_EXPORT static void SpreadCoarseToFine(const DataArrayDouble *coarseDA, const std::vector& coarseSt, DataArrayDouble *fineDA, const std::vector< std::pair >& fineLocInCoarse, const std::vector& facts); // diff --git a/src/MEDCoupling/MEDCouplingMemArray.hxx b/src/MEDCoupling/MEDCouplingMemArray.hxx index 38aa67d8f..7e01a1aaf 100644 --- a/src/MEDCoupling/MEDCouplingMemArray.hxx +++ b/src/MEDCoupling/MEDCouplingMemArray.hxx @@ -784,6 +784,7 @@ namespace ParaMEDMEM MEDCOUPLING_EXPORT void reprQuickOverview(std::ostream& stream) const; MEDCOUPLING_EXPORT void reprQuickOverviewData(std::ostream& stream, std::size_t maxNbOfByteInRepr) const; MEDCOUPLING_EXPORT bool isEqualIfNotWhy(const DataArrayChar& other, std::string& reason) const; + MEDCOUPLING_EXPORT std::vector toVectorOfBool() const; private: ~DataArrayByte() { } DataArrayByte() { } diff --git a/src/MEDCoupling/MEDCouplingMemArrayChar.cxx b/src/MEDCoupling/MEDCouplingMemArrayChar.cxx index 96e650be1..ce3baa868 100644 --- a/src/MEDCoupling/MEDCouplingMemArrayChar.cxx +++ b/src/MEDCoupling/MEDCouplingMemArrayChar.cxx @@ -2169,6 +2169,25 @@ bool DataArrayByte::isEqualIfNotWhy(const DataArrayChar& other, std::string& rea return DataArrayChar::isEqualIfNotWhy(other,reason); } +/*! + * This method is \b NOT wrapped into python because it can be useful only for performance reasons in C++ context. + * \throw if \a this is not allocated. + * \throw if \a this has not exactly one component. + */ +std::vector DataArrayByte::toVectorOfBool() const +{ + checkAllocated(); + if(getNumberOfComponents()!=1) + throw INTERP_KERNEL::Exception("DataArrayByte::toVectorOfBool : this method can be used only if this has one component !"); + int nbt(getNumberOfTuples()); + std::vector ret(nbt,false); + const char *pt(begin()); + for(int i=0;i >& partCompactFormat, int& axisId, int& sizeOfRange) +{ + int dim((int)partCompactFormat.size()); + int ret(-1); + for(int i=0;iret) + { + axisId=i; sizeOfRange=curDelta; + ret=curDelta; + } + } +} + +/*! + * This method is \b NOT wrapped in python because it has no sense in python (for performance reasons). + * This method starts from a structured mesh with structure \a st on which a boolean field \a crit is set. + * This method find for such minimalist information of mesh and field which is the part of the mesh, given by the range per axis in output parameter + * \a partCompactFormat that contains all the True in \a crit. The returned vector of boolean is the field reduced to that part. + * So the number of True is equal in \a st and in returned vector of boolean. + * + * \param [in] st - The structure per axis of the structured mesh considered. + * \param [in] crit - The field of boolean (for performance reasons) lying on the mesh defined by \a st. + * \param [out] partCompactFormat - The minimal part of \a st containing all the true of \a crit. + * \param [out] reducedCrit - The reduction of \a criterion on \a partCompactFormat. + * \return - The number of True in \a st (that is equal to those in \a reducedCrit) + */ +int MEDCouplingStructuredMesh::FindMinimalPartOf(const std::vector& st, const std::vector& crit, std::vector& reducedCrit, std::vector< std::pair >& partCompactFormat) +{ + if((int)crit.size()!=DeduceNumberOfGivenStructure(st)) + throw INTERP_KERNEL::Exception("MEDCouplingStructuredMesh::FindMinimalPartOf : size of vector of boolean is invalid regarding the declared structure !"); + switch((int)st.size()) + { + case 1: + { + return FindMinimalPartOf1D(st,crit,reducedCrit,partCompactFormat); + break; + } + case 2: + { + return FindMinimalPartOf2D(st,crit,reducedCrit,partCompactFormat); + break; + } + case 3: + { + return FindMinimalPartOf3D(st,crit,reducedCrit,partCompactFormat); + break; + } + default: + throw INTERP_KERNEL::Exception("MEDCouplingStructuredMesh::FindMinimalPartOf : only dimension 1, 2 and 3 are supported actually !"); + } +} + DataArrayInt *MEDCouplingStructuredMesh::Build1GTNodalConnectivity1D(const int *nodeStBg) { int nbOfCells(*nodeStBg-1); @@ -787,6 +846,134 @@ DataArrayInt *MEDCouplingStructuredMesh::Build1GTNodalConnectivityOfSubLevelMesh return conn.retn(); } +/*! + * \sa MEDCouplingStructuredMesh::FindMinimalPartOf + */ +int MEDCouplingStructuredMesh::FindMinimalPartOf1D(const std::vector& st, const std::vector& crit, std::vector& reducedCrit, std::vector< std::pair >& partCompactFormat) +{ + if(st.size()!=1) + throw INTERP_KERNEL::Exception("MEDCouplingStructuredMesh::FindMinimalPartOf1D : the input size of st must be equal to 1 !"); + int nxMin(std::numeric_limits::max()),nxMax(-std::numeric_limits::max()); + int nx(st[0]),ret(0); + for(int i=0;i& st, const std::vector& crit, std::vector& reducedCrit, std::vector< std::pair >& partCompactFormat) +{ + if(st.size()!=2) + throw INTERP_KERNEL::Exception("MEDCouplingStructuredMesh::FindMinimalPartOf2D : the input size of st must be equal to 2 !"); + int nxMin(std::numeric_limits::max()),nxMax(-std::numeric_limits::max()),nyMin(std::numeric_limits::max()),nyMax(-std::numeric_limits::max()); + int it(0),nx(st[0]),ny(st[1]); + int ret(0); + for(int i=0;i& st, const std::vector& crit, std::vector& reducedCrit, std::vector< std::pair >& partCompactFormat) +{ + if(st.size()!=3) + throw INTERP_KERNEL::Exception("MEDCouplingStructuredMesh::FindMinimalPartOf3D : the input size of st must be equal to 3 !"); + int nxMin(std::numeric_limits::max()),nxMax(-std::numeric_limits::max()),nyMin(std::numeric_limits::max()),nyMax(-std::numeric_limits::max()),nzMin(std::numeric_limits::max()),nzMax(-std::numeric_limits::max()); + int it(0),nx(st[0]),ny(st[1]),nz(st[2]); + int ret(0); + for(int i=0;i& st, const std::vector& crit, const std::vector< std::pair >& partCompactFormat, std::vector& reducedCrit) +{ + int nbt(DeduceNumberOfGivenRangeInCompactFrmt(partCompactFormat)); + std::vector dims(GetDimensionsFromCompactFrmt(partCompactFormat)); + reducedCrit.resize(nbt); + switch((int)st.size()) + { + case 1: + { + int nx(dims[0]); + int kk(partCompactFormat[0].first); + for(int i=0;i >& partCompactFormat); MEDCOUPLING_EXPORT static int DeduceNumberOfGivenStructure(const std::vector& st); + MEDCOUPLING_EXPORT static void FindTheWidestAxisOfGivenRangeInCompactFrmt(const std::vector< std::pair >& partCompactFormat, int& axisId, int& sizeOfRange); + MEDCOUPLING_EXPORT static int FindMinimalPartOf(const std::vector& st, const std::vector& crit, std::vector& reducedCrit, std::vector< std::pair >& partCompactFormat); private: static int GetNumberOfCellsOfSubLevelMesh(const std::vector& cgs, int mdim); static void GetReverseNodalConnectivity1(const std::vector& ngs, DataArrayInt *revNodal, DataArrayInt *revNodalIndx); @@ -89,6 +91,10 @@ namespace ParaMEDMEM static DataArrayInt *Build1GTNodalConnectivity3D(const int *nodeStBg); static DataArrayInt *Build1GTNodalConnectivityOfSubLevelMesh2D(const int *nodeStBg); static DataArrayInt *Build1GTNodalConnectivityOfSubLevelMesh3D(const int *nodeStBg); + static int FindMinimalPartOf1D(const std::vector& st, const std::vector& crit, std::vector& reducedCrit, std::vector< std::pair >& partCompactFormat); + static int FindMinimalPartOf2D(const std::vector& st, const std::vector& crit, std::vector& reducedCrit, std::vector< std::pair >& partCompactFormat); + static int FindMinimalPartOf3D(const std::vector& st, const std::vector& crit, std::vector& reducedCrit, std::vector< std::pair >& partCompactFormat); + static void ExtractVecOfBool(const std::vector& st, const std::vector& crit, const std::vector< std::pair >& partCompactFormat, std::vector& reducedCrit); protected: static int ZipNodeStructure(const int *nodeStBg, const int *nodeStEnd, int zipNodeSt[3]); protected: diff --git a/src/MEDCoupling_Swig/MEDCouplingBasicsTest.py b/src/MEDCoupling_Swig/MEDCouplingBasicsTest.py index de423c70c..c9e409eda 100644 --- a/src/MEDCoupling_Swig/MEDCouplingBasicsTest.py +++ b/src/MEDCoupling_Swig/MEDCouplingBasicsTest.py @@ -14861,7 +14861,7 @@ class MEDCouplingBasicsTest(unittest.TestCase): self.assertTrue(m2bis.isEqual(m2,1e-12)) # self.assertEqual(6,m2bis.getNumberOfCells())#3,2,4 - m2bis.refineWithFactor(3) + m2bis.refineWithFactor([3,3,3]) self.assertEqual(162,m2bis.getNumberOfCells()) self.assertEqual((7,4,10),m2bis.getNodeStruct()) self.assertEqual((1.5,3.5,2.5),m2bis.getOrigin()) @@ -14946,14 +14946,14 @@ class MEDCouplingBasicsTest(unittest.TestCase): self.assertEqual(0,amr.getNumberOfPatches()) self.assertEqual(1,amr.getMaxNumberOfLevelsRelativeToThis()) self.assertEqual(2,amr.getSpaceDimension()) - amr.addPatch([(1,2),(0,1)],4) + amr.addPatch([(1,2),(0,1)],[4,4]) self.assertEqual(4,amr.getNumberOfCellsAtCurrentLevel()) self.assertEqual(20,amr.getNumberOfCellsRecursiveWithOverlap()) self.assertEqual(19,amr.getNumberOfCellsRecursiveWithoutOverlap()) self.assertEqual(1,amr.getNumberOfPatches()) self.assertEqual(2,amr.getMaxNumberOfLevelsRelativeToThis()) self.assertEqual(2,amr.getSpaceDimension()) - amr[0].addPatch([(2,3),(1,3)],2) + amr[0].addPatch([(2,3),(1,3)],[2,2]) self.assertEqual(amr[0].getBLTRRange(),[(1,2),(0,1)]) self.assertEqual(4,amr.getNumberOfCellsAtCurrentLevel()) self.assertEqual(28,amr.getNumberOfCellsRecursiveWithOverlap()) @@ -14961,7 +14961,7 @@ class MEDCouplingBasicsTest(unittest.TestCase): self.assertEqual(1,amr.getNumberOfPatches()) self.assertEqual(3,amr.getMaxNumberOfLevelsRelativeToThis()) self.assertEqual(2,amr.getSpaceDimension()) - amr[0].addPatch([(0,2),(3,4)],3) + amr[0].addPatch([(0,2),(3,4)],[3,3]) self.assertEqual(16,amr[0].getMesh().getNumberOfCellsAtCurrentLevel()) self.assertEqual(46,amr.getNumberOfCellsRecursiveWithOverlap()) self.assertEqual(41,amr.getNumberOfCellsRecursiveWithoutOverlap()) diff --git a/src/MEDCoupling_Swig/MEDCouplingCommon.i b/src/MEDCoupling_Swig/MEDCouplingCommon.i index 1f8160779..4c1ba6391 100644 --- a/src/MEDCoupling_Swig/MEDCouplingCommon.i +++ b/src/MEDCoupling_Swig/MEDCouplingCommon.i @@ -45,6 +45,7 @@ #include "MEDCouplingTypemaps.i" #include "InterpKernelAutoPtr.hxx" +#include "BoxSplittingOptions.hxx" using namespace ParaMEDMEM; using namespace INTERP_KERNEL; @@ -367,6 +368,30 @@ using namespace INTERP_KERNEL; %include "MEDCouplingRefCountObject.i" %include "MEDCouplingMemArray.i" +namespace INTERP_KERNEL +{ + /*! + * \class BoxSplittingOptions + * Class defining the options for box splitting used for AMR algorithm like creation of patches following a criterion. + */ + class BoxSplittingOptions + { + public: + BoxSplittingOptions(); + void init() throw(INTERP_KERNEL::Exception); + double getEffeciency() const throw(INTERP_KERNEL::Exception); + void setEffeciency(double effeciency) throw(INTERP_KERNEL::Exception); + double getEffeciencySnd() const throw(INTERP_KERNEL::Exception); + void setEffeciencySnd(double effeciencySnd) throw(INTERP_KERNEL::Exception); + int getMinCellDirection() const throw(INTERP_KERNEL::Exception); + void setMinCellDirection(int minCellDirection) throw(INTERP_KERNEL::Exception); + int getMaxCells() const throw(INTERP_KERNEL::Exception); + void setMaxCells(int maxCells) throw(INTERP_KERNEL::Exception); + void copyOptions(const BoxSplittingOptions & other) throw(INTERP_KERNEL::Exception); + std::string printOptions() const throw(INTERP_KERNEL::Exception); + }; +} + namespace ParaMEDMEM { typedef enum @@ -2994,7 +3019,7 @@ namespace ParaMEDMEM std::string getAxisUnit() const throw(INTERP_KERNEL::Exception); double getMeasureOfAnyCell() const throw(INTERP_KERNEL::Exception); MEDCouplingCMesh *convertToCartesian() const throw(INTERP_KERNEL::Exception); - void refineWithFactor(int factor) throw(INTERP_KERNEL::Exception); + void refineWithFactor(const std::vector& factors) throw(INTERP_KERNEL::Exception); %extend { MEDCouplingIMesh() @@ -4632,11 +4657,11 @@ namespace ParaMEDMEM return ret; } - void addPatch(PyObject *bottomLeftTopRight, int factor) throw(INTERP_KERNEL::Exception) + void addPatch(PyObject *bottomLeftTopRight, const std::vector& factors) throw(INTERP_KERNEL::Exception) { std::vector< std::pair > inp; convertPyToVectorPairInt(bottomLeftTopRight,inp); - self->addPatch(inp,factor); + self->addPatch(inp,factors); } MEDCouplingCartesianAMRPatch *__getitem__(int patchId) const throw(INTERP_KERNEL::Exception) @@ -4715,11 +4740,11 @@ namespace ParaMEDMEM return ParaMEDMEM_MEDCouplingCartesianAMRMesh_New(meshName,spaceDim,nodeStrct,origin,dxyz); } - void addPatch(PyObject *bottomLeftTopRight, int factor) throw(INTERP_KERNEL::Exception) + void addPatch(PyObject *bottomLeftTopRight, const std::vector& factors) throw(INTERP_KERNEL::Exception) { std::vector< std::pair > inp; convertPyToVectorPairInt(bottomLeftTopRight,inp); - self->addPatch(inp,factor); + self->addPatch(inp,factors); } MEDCouplingCartesianAMRMesh *getFather() const throw(INTERP_KERNEL::Exception) diff --git a/src/MEDCoupling_Swig/MEDCouplingMemArray.i b/src/MEDCoupling_Swig/MEDCouplingMemArray.i index e7e283efc..27f0d46cf 100644 --- a/src/MEDCoupling_Swig/MEDCouplingMemArray.i +++ b/src/MEDCoupling_Swig/MEDCouplingMemArray.i @@ -5194,6 +5194,223 @@ namespace ParaMEDMEM return ParaMEDMEM_DataArrayByte_presenceOfTuple(self,obj); } } + + DataArrayByte *__setitem__(PyObject *obj, PyObject *value) throw(INTERP_KERNEL::Exception) + { + self->checkAllocated(); + const char msg[]="Unexpected situation in __setitem__ !"; + int nbOfTuples(self->getNumberOfTuples()),nbOfComponents(self->getNumberOfComponents()); + int sw1,sw2; + int i1; + std::vector v1; + DataArrayInt *d1=0; + DataArrayIntTuple *dd1=0; + convertObjToPossibleCpp1(value,sw1,i1,v1,d1,dd1); + int it1,ic1; + std::vector vt1,vc1; + std::pair > pt1,pc1; + DataArrayInt *dt1=0,*dc1=0; + convertObjToPossibleCpp3(obj,nbOfTuples,nbOfComponents,sw2,it1,ic1,vt1,vc1,pt1,pc1,dt1,dc1); + MEDCouplingAutoRefCountObjectPtr tmp; + switch(sw2) + { + case 1: + { + switch(sw1) + { + case 1: + self->setPartOfValuesSimple1(i1,it1,it1+1,1,0,nbOfComponents,1); + return self; + default: + throw INTERP_KERNEL::Exception(msg); + } + break; + } + case 2: + { + switch(sw1) + { + case 1: + self->setPartOfValuesSimple3(i1,&vt1[0],&vt1[0]+vt1.size(),0,nbOfComponents,1); + return self; + default: + throw INTERP_KERNEL::Exception(msg); + } + break; + } + case 3: + { + switch(sw1) + { + case 1: + self->setPartOfValuesSimple1(i1,pt1.first,pt1.second.first,pt1.second.second,0,nbOfComponents,1); + return self; + default: + throw INTERP_KERNEL::Exception(msg); + } + break; + } + case 4: + { + switch(sw1) + { + case 1: + self->setPartOfValuesSimple3(i1,dt1->getConstPointer(),dt1->getConstPointer()+dt1->getNbOfElems(),0,nbOfComponents,1); + return self; + default: + throw INTERP_KERNEL::Exception(msg); + } + break; + } + case 5: + { + switch(sw1) + { + case 1: + self->setPartOfValuesSimple1(i1,it1,it1+1,1,ic1,ic1+1,1); + return self; + default: + throw INTERP_KERNEL::Exception(msg); + } + break; + } + case 6: + { + switch(sw1) + { + case 1: + self->setPartOfValuesSimple3(i1,&vt1[0],&vt1[0]+vt1.size(),ic1,ic1+1,1); + return self; + default: + throw INTERP_KERNEL::Exception(msg); + } + break; + } + case 7: + { + switch(sw1) + { + case 1: + self->setPartOfValuesSimple1(i1,pt1.first,pt1.second.first,pt1.second.second,ic1,ic1+1,1); + return self; + default: + throw INTERP_KERNEL::Exception(msg); + } + break; + } + case 8: + { + switch(sw1) + { + case 1: + self->setPartOfValuesSimple3(i1,dt1->getConstPointer(),dt1->getConstPointer()+dt1->getNbOfElems(),ic1,ic1+1,1); + return self; + default: + throw INTERP_KERNEL::Exception(msg); + } + break; + } + case 9: + { + switch(sw1) + { + case 1: + self->setPartOfValuesSimple2(i1,&it1,&it1+1,&vc1[0],&vc1[0]+vc1.size()); + return self; + default: + throw INTERP_KERNEL::Exception(msg); + } + break; + } + case 10: + { + switch(sw1) + { + case 1: + self->setPartOfValuesSimple2(i1,&vt1[0],&vt1[0]+vt1.size(),&vc1[0],&vc1[0]+vc1.size()); + return self; + default: + throw INTERP_KERNEL::Exception(msg); + } + break; + } + case 11: + { + switch(sw1) + { + case 1: + self->setPartOfValuesSimple4(i1,pt1.first,pt1.second.first,pt1.second.second,&vc1[0],&vc1[0]+vc1.size()); + return self; + default: + throw INTERP_KERNEL::Exception(msg); + } + break; + } + case 12: + { + switch(sw1) + { + case 1: + self->setPartOfValuesSimple2(i1,dt1->getConstPointer(),dt1->getConstPointer()+dt1->getNbOfElems(),&vc1[0],&vc1[0]+vc1.size()); + return self; + default: + throw INTERP_KERNEL::Exception(msg); + } + break; + } + case 13: + { + switch(sw1) + { + case 1: + self->setPartOfValuesSimple1(i1,it1,it1+1,1,pc1.first,pc1.second.first,pc1.second.second); + return self; + default: + throw INTERP_KERNEL::Exception(msg); + } + break; + } + case 14: + { + switch(sw1) + { + case 1: + self->setPartOfValuesSimple3(i1,&vt1[0],&vt1[0]+vt1.size(),pc1.first,pc1.second.first,pc1.second.second); + return self; + default: + throw INTERP_KERNEL::Exception(msg); + } + break; + } + case 15: + { + switch(sw1) + { + case 1: + self->setPartOfValuesSimple1(i1,pt1.first,pt1.second.first,pt1.second.second,pc1.first,pc1.second.first,pc1.second.second); + return self; + default: + throw INTERP_KERNEL::Exception(msg); + } + break; + } + case 16: + { + switch(sw1) + { + case 1: + self->setPartOfValuesSimple3(i1,dt1->getConstPointer(),dt1->getConstPointer()+dt1->getNbOfElems(),pc1.first,pc1.second.first,pc1.second.second); + return self; + default: + throw INTERP_KERNEL::Exception(msg); + } + break; + } + default: + throw INTERP_KERNEL::Exception(msg); + } + return self; + } } };