From: geay Date: Wed, 21 May 2014 15:42:23 +0000 (+0200) Subject: First test OK with createPatchesFromCriterion X-Git-Tag: V7_5_0a1~37^2~40 X-Git-Url: http://git.salome-platform.org/gitweb/?a=commitdiff_plain;h=05c2a9751ce5dbf5db11fa2153f802416f8ddf50;p=modules%2Fmed.git First test OK with createPatchesFromCriterion --- diff --git a/src/MEDCoupling/MEDCouplingCartesianAMRMesh.cxx b/src/MEDCoupling/MEDCouplingCartesianAMRMesh.cxx index 929adf7d6..d97d294cf 100644 --- a/src/MEDCoupling/MEDCouplingCartesianAMRMesh.cxx +++ b/src/MEDCoupling/MEDCouplingCartesianAMRMesh.cxx @@ -19,11 +19,13 @@ // Author : Anthony Geay #include "MEDCouplingCartesianAMRMesh.hxx" +#include "MEDCoupling1GTUMesh.hxx" #include "MEDCouplingIMesh.hxx" #include "MEDCouplingUMesh.hxx" #include #include +#include using namespace ParaMEDMEM; @@ -179,81 +181,270 @@ public: int getNumberOfCells() const { return (int)_crit.size(); } void setNumberOfTrue(int nboft) { _nb_of_true=nboft; } std::vector& getCriterion() { return _crit; } + const std::vector& getConstCriterion() const { return _crit; } + void setPart(const std::vector< std::pair >& part) { _part=part; } std::vector< std::pair >& getPart() { return _part; } const std::vector< std::pair >& getConstPart() const { return _part; } bool presenceOfTrue() const { return _nb_of_true>0; } + std::vector computeCGS() const { return MEDCouplingStructuredMesh::GetDimensionsFromCompactFrmt(_part); } + std::vector< std::vector > computeSignature() const { return MEDCouplingStructuredMesh::ComputeSignaturePerAxisOf(computeCGS(),getConstCriterion()); } + double getEfficiencyPerAxis(int axisId) const { return (double)_nb_of_true/((double)(_part[axisId].second-_part[axisId].first)); } + void zipToFitOnCriterion(); + void updateNumberOfTrue() const; + MEDCouplingAutoRefCountObjectPtr extractPart(const std::vector< std::pair >&partInGlobal) const; + MEDCouplingAutoRefCountObjectPtr deepCpy() const; protected: ~InternalPatch() { } private: - int _nb_of_true; + mutable int _nb_of_true; std::vector _crit; + //! _part is global 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 +void InternalPatch::zipToFitOnCriterion() { - int cut_found = 0 ; - int cut_place = -1 ; - float * ratio = NULL ; - float * ratio_inner = NULL ; + std::vector cgs(computeCGS()); + std::vector newCrit; + std::vector< std::pair > newPart,newPart2; + int newNbOfTrue(MEDCouplingStructuredMesh::FindMinimalPartOf(cgs,_crit,newCrit,newPart)); + MEDCouplingStructuredMesh::ChangeReferenceToGlobalOfCompactFrmt(_part,newPart,newPart2); + if(newNbOfTrue!=_nb_of_true) + throw INTERP_KERNEL::Exception("InternalPatch::zipToFitOnCrit : internal error !"); + _crit=newCrit; _part=newPart2; +} + +void InternalPatch::updateNumberOfTrue() const +{ + _nb_of_true=(int)std::count(_crit.begin(),_crit.end(),true); +} + +MEDCouplingAutoRefCountObjectPtr InternalPatch::extractPart(const std::vector< std::pair >&partInGlobal) const +{ + MEDCouplingAutoRefCountObjectPtr ret(new InternalPatch); + std::vector cgs(computeCGS()); + std::vector< std::pair > newPart; + MEDCouplingStructuredMesh::ChangeReferenceFromGlobalOfCompactFrmt(_part,partInGlobal,newPart); + MEDCouplingStructuredMesh::ExtractFieldOfBoolFrom(cgs,_crit,newPart,ret->getCriterion()); + ret->setPart(partInGlobal); + ret->updateNumberOfTrue(); + return ret; +} + +MEDCouplingAutoRefCountObjectPtr InternalPatch::deepCpy() const +{ + MEDCouplingAutoRefCountObjectPtr ret(new InternalPatch); + (*ret)=*this; + return ret; +} - ratio = new float [big_dims-1]; - for(unsigned int id=0; id ratio(rangeOfAxisId-1); + for(int id=0;id > rectH(patchToBeSplit->getConstPart()); + if(h==0) + rectH[axisId].second=patchToBeSplit->getConstPart()[axisId].first+id; else - nb_cells_h = patch_h.getNy() ; + rectH[axisId].first=patchToBeSplit->getConstPart()[axisId].first+id; + MEDCouplingAutoRefCountObjectPtr p(patchToBeSplit->deepCpy()); + p->zipToFitOnCriterion(); + //anouar rectH ? + efficiency[h]=p->getEfficiencyPerAxis(axisId); + } + ratio[id]=std::max(efficiency[0],efficiency[1])/std::min(efficiency[0],efficiency[1]); + } + int minCellDirection(bso.getMinCellDirection()),indexMin(-1); + int dimRatioInner(rangeOfAxisId-1-2*(minCellDirection-1)); + std::vector ratio_inner(dimRatioInner); + double minRatio(1.e10); + for(int i=0; igetConstPart()[axisId].first-1; +} - int nb_cells_flag_h = patch_h.getNumberOfCellsFlags(); - efficiency[h] = float (nb_cells_flag_h) / float(nb_cells_h) ; +void FindHole(const INTERP_KERNEL::BoxSplittingOptions& bso, const InternalPatch *patchToBeSplit, const int axisId, bool& cutFound, int& cutPlace) +{ + cutPlace=-1; cutFound=false; + int minCellDirection(bso.getMinCellDirection()); + + int sortedDims[2]; + sortedDims[0]=axisId==0?1:0; + sortedDims[1]=axisId==0?0:1; + static const int dim(2); + std::vector< std::vector > signatures(patchToBeSplit->computeSignature()); + for(int id=0;id& signature(signatures[sortedDims[id]]); + std::vector hole; + std::vector distance ; + int len((int)signature.size()); + for(int i=0;i= 2*minCellDirection && i >= minCellDirection-1 && i <= len-minCellDirection-1) + hole.push_back(i); + if (hole.size()>0) + { + double center(((double)len/2.)+0.5); + for(std::size_t i=0;igetConstPart()[axisId].first+1; + return ; } - ratio[id] = max(efficiency[0],efficiency[1])/ - min(efficiency[0],efficiency[1]) ; } +} + +void FindInflection(const INTERP_KERNEL::BoxSplittingOptions& bso, const InternalPatch *patchToBeSplit, bool& cutFound, int& cutPlace, int& axisId) +{ + cutFound=false; cutPlace=-1; axisId=-1; - int dim_ratio_inner = big_dims-1-2*(minCellDirection-1) ; - ratio_inner = new float [dim_ratio_inner]; - float min_ratio = 1.E10 ; - int index_min = -1 ; - for(int i=0; i >& part(patchToBeSplit->getConstPart()); + int sign,minCellDirection(bso.getMinCellDirection()); + static const int dim = 2 ; + + std::vector zeroCrossDims(2,-1); + std::vector zeroCrossVals(2,-1); + std::vector< std::vector > signatures(patchToBeSplit->computeSignature()); + for (int id=0;id& signature(signatures[id]); + + std::vector derivate_second_order,gradient_absolute,signe_change,zero_cross,edge,max_cross_list ; + std::vector distance ; + + for (unsigned int i=1;i 0 ) + sign = 1 ; + if (derivate_second_order[i]*derivate_second_order[i+1] == 0 ) + sign = 0 ; + if ( sign==0 || sign==-1 ) + if ( i >= (unsigned int)minCellDirection-2 && i <= signature.size()-minCellDirection-2 ) + { + zero_cross.push_back(i) ; + edge.push_back(gradient_absolute[i]) ; + } + signe_change.push_back(sign) ; } + if ( zero_cross.size() > 0 ) + { + int max_cross=*max_element(edge.begin(),edge.end()) ; + for (unsigned int i=0;i=0 ) + { + zeroCrossDims[id] = best_place ; + zeroCrossVals[id] = max_cross ; + } + } + derivate_second_order.clear() ; + gradient_absolute.clear() ; + signe_change.clear() ; + zero_cross.clear() ; + edge.clear() ; + max_cross_list.clear() ; + distance.clear() ; } - cut_found = 1 ; - cut_place = index_min + _indexCorners[dissect_direction] - 1 ; - cut[0] = cut_found ; - cut[1] = cut_place ; - cut[2] = dissect_direction ; - delete [] ratio ; - delete [] ratio_inner ; -} -#endif + + if ( zeroCrossDims[0]!=-1 || zeroCrossDims[1]!=-1 ) + { + int max_cross_dims = *max_element(zeroCrossVals.begin(),zeroCrossVals.end()) ; + + if (zeroCrossVals[0]==max_cross_dims && zeroCrossVals[1]==max_cross_dims ) + { + int nl_left(part[0].second-part[0].first); + int nc_left(part[1].second-part[1].first); + if ( nl_left >= nc_left ) + max_cross_dims = 0 ; + else + max_cross_dims = 1 ; + } + else + max_cross_dims=std::find(zeroCrossVals.begin(),zeroCrossVals.end(),max_cross_dims)-zeroCrossVals.begin(); + cutFound=true; + cutPlace=zeroCrossDims[max_cross_dims]; + axisId=max_cross_dims ; + } +} + +void TryAction4(const INTERP_KERNEL::BoxSplittingOptions& bso, const InternalPatch *patchToBeSplit, int axisId, int rangeOfAxisId, bool& cutFound, int& cutPlace) +{ + cutFound=false; + if(patchToBeSplit->getEfficiency() <= bso.getEffeciencySnd()) + { + if(rangeOfAxisId>=2*bso.getMinCellDirection()) + { + cutFound=true; + cutPlace=rangeOfAxisId/2+patchToBeSplit->getConstPart()[axisId].first-1; + } + } + else + { + if(patchToBeSplit->getNumberOfCells()>bso.getMaxCells()) + { + DissectBigPatch(bso,patchToBeSplit,axisId,rangeOfAxisId,cutFound,cutPlace); + } + } +} + +MEDCouplingAutoRefCountObjectPtr DealWithNoCut(const InternalPatch *patch) +{ + MEDCouplingAutoRefCountObjectPtr ret(const_cast(patch)); + ret->incrRef(); + return ret; +} + +void DealWithCut(const InternalPatch *patchToBeSplit, int axisId, int cutPlace, std::vector >& listOfPatches) +{ + MEDCouplingAutoRefCountObjectPtr leftPart,rightPart; + std::vector< std::pair > rect(patchToBeSplit->getConstPart()); + std::vector< std::pair > leftRect(rect),rightRect(rect); + leftRect[axisId].second=cutPlace+1; + rightRect[axisId].first=cutPlace+1; + leftPart=patchToBeSplit->extractPart(leftRect); + rightPart=patchToBeSplit->extractPart(rightRect); + leftPart->zipToFitOnCriterion(); rightPart->zipToFitOnCriterion(); + listOfPatches.push_back(leftPart); + listOfPatches.push_back(rightPart); +} + /// @endcond /*! - * This method creates patches in \a this (by destroying the patches if any). This method uses \a criterion + * This method creates patches in \a this (by destroying the patches if any). This method uses \a criterion array as a field on cells on this level. */ void MEDCouplingCartesianAMRMesh::createPatchesFromCriterion(const INTERP_KERNEL::BoxSplittingOptions& bso, const DataArrayByte *criterion, const std::vector& factors) { @@ -278,13 +469,22 @@ void MEDCouplingCartesianAMRMesh::createPatchesFromCriterion(const INTERP_KERNEL for(std::vector< MEDCouplingAutoRefCountObjectPtr >::iterator it=listOfPatches.begin();it!=listOfPatches.end();it++) { // - if((*it)->getEfficiency()>=bso.getEffeciency()) - { - if((*it)->getNumberOfCells()>=bso.getMaxCells()) - { - - } - } + int axisId,rangeOfAxisId; + bool cutFound; + int cutPlace; + MEDCouplingStructuredMesh::FindTheWidestAxisOfGivenRangeInCompactFrmt((*it)->getConstPart(),axisId,rangeOfAxisId); + if((*it)->getEfficiency()>=bso.getEffeciency() && (*it)->getNumberOfCells() cells; + std::vector< MEDCouplingAutoRefCountObjectPtr > cellsSafe; + for(std::vector< MEDCouplingAutoRefCountObjectPtr >::const_iterator it=_patches.begin();it!=_patches.end();it++) + { + const MEDCouplingCartesianAMRPatch *patch(*it); + if(patch) + { + MEDCouplingAutoRefCountObjectPtr cell(patch->getMesh()->getImageMesh()->asSingleCell()); + MEDCouplingAutoRefCountObjectPtr cell1SGT(cell->build1SGTUnstructured()); + cellsSafe.push_back(cell1SGT); cells.push_back(cell1SGT); + } + } + return MEDCoupling1SGTUMesh::Merge1SGTUMeshes(cells); +} + MEDCouplingCartesianAMRMesh::MEDCouplingCartesianAMRMesh(const std::string& meshName, int spaceDim, const int *nodeStrctStart, const int *nodeStrctStop, const double *originStart, const double *originStop, const double *dxyzStart, const double *dxyzStop):_father(0) { diff --git a/src/MEDCoupling/MEDCouplingCartesianAMRMesh.hxx b/src/MEDCoupling/MEDCouplingCartesianAMRMesh.hxx index 5f6309f48..d0ab4a421 100644 --- a/src/MEDCoupling/MEDCouplingCartesianAMRMesh.hxx +++ b/src/MEDCoupling/MEDCouplingCartesianAMRMesh.hxx @@ -34,6 +34,7 @@ namespace ParaMEDMEM class MEDCouplingIMesh; class MEDCouplingUMesh; class DataArrayByte; + class MEDCoupling1SGTUMesh; class MEDCouplingCartesianAMRMesh; /// @cond INTERNAL @@ -77,6 +78,7 @@ namespace ParaMEDMEM MEDCOUPLING_EXPORT int getNumberOfCellsAtCurrentLevel() const; MEDCOUPLING_EXPORT int getNumberOfCellsRecursiveWithOverlap() const; MEDCOUPLING_EXPORT int getNumberOfCellsRecursiveWithoutOverlap() const; + MEDCOUPLING_EXPORT const MEDCouplingIMesh *getImageMesh() const { return _mesh; } // MEDCOUPLING_EXPORT const MEDCouplingCartesianAMRMesh *getFather() const; MEDCOUPLING_EXPORT const MEDCouplingCartesianAMRMesh *getGodFather() const; @@ -88,6 +90,7 @@ namespace ParaMEDMEM MEDCOUPLING_EXPORT const MEDCouplingCartesianAMRPatch *getPatch(int patchId) const; // MEDCOUPLING_EXPORT MEDCouplingUMesh *buildUnstructured() const; + MEDCOUPLING_EXPORT MEDCoupling1SGTUMesh *buildMeshFromPatchEnvelop() const; private: MEDCouplingCartesianAMRMesh(const std::string& meshName, int spaceDim, const int *nodeStrctStart, const int *nodeStrctStop, const double *originStart, const double *originStop, const double *dxyzStart, const double *dxyzStop); diff --git a/src/MEDCoupling/MEDCouplingIMesh.cxx b/src/MEDCoupling/MEDCouplingIMesh.cxx index 142fac4b8..4b3da324a 100644 --- a/src/MEDCoupling/MEDCouplingIMesh.cxx +++ b/src/MEDCoupling/MEDCouplingIMesh.cxx @@ -206,6 +206,36 @@ void MEDCouplingIMesh::refineWithFactor(const std::vector& factors) declareAsNew(); } +/*! + * This method returns a newly created mesh containing a single cell in it. This returned cell covers exactly the space covered by \a this. + * + * \return MEDCouplingIMesh * - A newly created object (to be managed by the caller with decrRef) containing simply one cell. + * + * \throw if \a this does not pass the \c checkCoherency test. + */ +MEDCouplingIMesh *MEDCouplingIMesh::asSingleCell() const +{ + checkCoherency(); + int spaceDim(getSpaceDimension()),nodeSt[3]; + double dxyz[3]; + for(int i=0;i=2) + { + nodeSt[i]=2; + dxyz[i]=(_structure[i]-1)*_dxyz[i]; + } + else + { + nodeSt[i]=_structure[i]; + dxyz[i]=_dxyz[i]; + } + } + MEDCouplingAutoRefCountObjectPtr ret(MEDCouplingIMesh::New(getName(),getSpaceDimension(),nodeSt,nodeSt+spaceDim,_origin,_origin+spaceDim,dxyz,dxyz+spaceDim)); + ret->copyTinyInfoFrom(this); + return ret.retn(); +} + /*! * This static method is useful to condense field on cells of a MEDCouplingIMesh instance coming from a refinement ( MEDCouplingIMesh::refineWithFactor for example) * to a coarse MEDCouplingIMesh instance. So this method can be seen as a specialization in P0P0 conservative interpolation non overlaping from fine image mesh diff --git a/src/MEDCoupling/MEDCouplingIMesh.hxx b/src/MEDCoupling/MEDCouplingIMesh.hxx index 896a7835d..ad6141254 100644 --- a/src/MEDCoupling/MEDCouplingIMesh.hxx +++ b/src/MEDCoupling/MEDCouplingIMesh.hxx @@ -47,6 +47,7 @@ namespace ParaMEDMEM MEDCOUPLING_EXPORT double getMeasureOfAnyCell() const; MEDCOUPLING_EXPORT MEDCouplingCMesh *convertToCartesian() const; MEDCOUPLING_EXPORT void refineWithFactor(const std::vector& factors); + MEDCOUPLING_EXPORT MEDCouplingIMesh *asSingleCell() const; 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/MEDCouplingStructuredMesh.cxx b/src/MEDCoupling/MEDCouplingStructuredMesh.cxx index 842cba3ed..19c054b6e 100644 --- a/src/MEDCoupling/MEDCouplingStructuredMesh.cxx +++ b/src/MEDCoupling/MEDCouplingStructuredMesh.cxx @@ -739,26 +739,79 @@ int MEDCouplingStructuredMesh::FindMinimalPartOf(const std::vector& st, con { if((int)crit.size()!=DeduceNumberOfGivenStructure(st)) throw INTERP_KERNEL::Exception("MEDCouplingStructuredMesh::FindMinimalPartOf : size of vector of boolean is invalid regarding the declared structure !"); + int ret(-1); switch((int)st.size()) { case 1: { - return FindMinimalPartOf1D(st,crit,reducedCrit,partCompactFormat); + ret=FindMinimalPartOf1D(st,crit,partCompactFormat); break; } case 2: { - return FindMinimalPartOf2D(st,crit,reducedCrit,partCompactFormat); + ret=FindMinimalPartOf2D(st,crit,partCompactFormat); break; } case 3: { - return FindMinimalPartOf3D(st,crit,reducedCrit,partCompactFormat); + ret=FindMinimalPartOf3D(st,crit,partCompactFormat); break; } default: throw INTERP_KERNEL::Exception("MEDCouplingStructuredMesh::FindMinimalPartOf : only dimension 1, 2 and 3 are supported actually !"); } + ExtractFieldOfBoolFrom(st,crit,partCompactFormat,reducedCrit); + return ret; +} + +/*! + * This method is \b NOT wrapped in python. + * This method considers \a crit input parameter as a matrix having dimensions specified by \a st. This method returns for each axis + * the signature, that is to say the number of elems equal to true in \a crit along this axis. + */ +std::vector< std::vector > MEDCouplingStructuredMesh::ComputeSignaturePerAxisOf(const std::vector& st, const std::vector& crit) +{ + int dim((int)st.size()); + std::vector< std::vector > ret(dim); + switch(dim) + { + case 1: + { + int nx(st[0]); + ret[0].resize(nx); + std::vector& retX(ret[0]); + for(int i=0;i& retX(ret[0]); + for(int i=0;i& retY(ret[1]); + for(int j=0;j& st, const std::vector& crit, std::vector& reducedCrit, std::vector< std::pair >& partCompactFormat) +int MEDCouplingStructuredMesh::FindMinimalPartOf1D(const std::vector& st, const std::vector& crit, 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 !"); @@ -867,14 +920,13 @@ int MEDCouplingStructuredMesh::FindMinimalPartOf1D(const std::vector& st, c return ret; partCompactFormat.resize(1); partCompactFormat[0].first=nxMin; partCompactFormat[0].second=nxMax+1; - ExtractVecOfBool(st,crit,partCompactFormat,reducedCrit); return ret; } /*! * \sa MEDCouplingStructuredMesh::FindMinimalPartOf */ -int MEDCouplingStructuredMesh::FindMinimalPartOf2D(const std::vector& st, const std::vector& crit, std::vector& reducedCrit, std::vector< std::pair >& partCompactFormat) +int MEDCouplingStructuredMesh::FindMinimalPartOf2D(const std::vector& st, const std::vector& crit, 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 !"); @@ -896,14 +948,13 @@ int MEDCouplingStructuredMesh::FindMinimalPartOf2D(const std::vector& st, c partCompactFormat.resize(2); partCompactFormat[0].first=nxMin; partCompactFormat[0].second=nxMax+1; partCompactFormat[1].first=nyMin; partCompactFormat[1].second=nyMax+1; - ExtractVecOfBool(st,crit,partCompactFormat,reducedCrit); return ret; } /*! * \sa MEDCouplingStructuredMesh::FindMinimalPartOf */ -int MEDCouplingStructuredMesh::FindMinimalPartOf3D(const std::vector& st, const std::vector& crit, std::vector& reducedCrit, std::vector< std::pair >& partCompactFormat) +int MEDCouplingStructuredMesh::FindMinimalPartOf3D(const std::vector& st, const std::vector& crit, 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 !"); @@ -928,7 +979,6 @@ int MEDCouplingStructuredMesh::FindMinimalPartOf3D(const std::vector& st, c partCompactFormat[0].first=nxMin; partCompactFormat[0].second=nxMax+1; partCompactFormat[1].first=nyMin; partCompactFormat[1].second=nyMax+1; partCompactFormat[2].first=nzMin; partCompactFormat[2].second=nzMax+1; - ExtractVecOfBool(st,crit,partCompactFormat,reducedCrit); return ret; } @@ -1234,7 +1284,7 @@ std::vector MEDCouplingStructuredMesh::GetDimensionsFromCompactFrmt(const s * \param [in] partCompactFormat The compact subpart to be enabled. * \param [in,out] vectToSwitchOn Vector which fetched items are enabled. * - * \sa MEDCouplingStructuredMesh::BuildExplicitIdsFrom + * \sa MEDCouplingStructuredMesh::BuildExplicitIdsFrom, ExtractFieldOfBoolFrom */ void MEDCouplingStructuredMesh::SwitchOnIdsFrom(const std::vector& st, const std::vector< std::pair >& partCompactFormat, std::vector& vectToSwitchOn) { @@ -1276,16 +1326,144 @@ void MEDCouplingStructuredMesh::SwitchOnIdsFrom(const std::vector& st, cons break; } default: - throw INTERP_KERNEL::Exception("MEDCouplingStructuredMesh::BuildExplicitIdsFrom : Dimension supported are 1,2 or 3 !"); + throw INTERP_KERNEL::Exception("MEDCouplingStructuredMesh::SwitchOnIdsFrom : Dimension supported are 1,2 or 3 !"); + } +} + +/*! + * Obviously this method is \b NOT wrapped in python. + * This method is close to SwitchOnIdsFrom except that here, a sub field \a fieldOut is built starting from the input field \a fieldOfBool having the structure \a st. + * The extraction is defined by \a partCompactFormat. + * + * \param [in] st The entity structure. + * \param [in] fieldOfBool field of booleans having the size equal to \c MEDCouplingStructuredMesh::DeduceNumberOfGivenStructure(st). + * \param [in] partCompactFormat The compact subpart to be enabled. + * \param [out] fieldOut the result of the extraction. + * + * \sa MEDCouplingStructuredMesh::BuildExplicitIdsFrom, SwitchOnIdsFrom + */ +void MEDCouplingStructuredMesh::ExtractFieldOfBoolFrom(const std::vector& st, const std::vector& fieldOfBool, const std::vector< std::pair >& partCompactFormat, std::vector& fieldOut) +{ + if(st.size()!=partCompactFormat.size()) + throw INTERP_KERNEL::Exception("MEDCouplingStructuredMesh::ExtractFieldOfBoolFrom : input arrays must have the same size !"); + if((int)fieldOfBool.size()!=DeduceNumberOfGivenStructure(st)) + throw INTERP_KERNEL::Exception("MEDCouplingStructuredMesh::ExtractFieldOfBoolFrom : invalid size of input field of boolean regarding the structure !"); + std::vector dims(GetDimensionsFromCompactFrmt(partCompactFormat)); + int nbOfTuplesOfOutField(DeduceNumberOfGivenStructure(dims)); + fieldOut.resize(nbOfTuplesOfOutField); + int it(0); + switch(st.size()) + { + case 3: + { + for(int i=0;i >& bigInAbs, const std::vector< std::pair >& partOfBigInAbs, std::vector< std::pair >& partOfBigRelativeToBig) +{ + std::size_t dim(bigInAbs.size()); + if(dim!=partOfBigInAbs.size()) + throw INTERP_KERNEL::Exception("MEDCouplingStructuredMesh::ChangeReferenceFromGlobalOfCompactFrmt : The size of parts (dimension) must be the same !"); + partOfBigRelativeToBig.resize(dim); + for(std::size_t i=0;ibigInAbs[i].second) + { + std::ostringstream oss; oss << "MEDCouplingStructuredMesh::ChangeReferenceFromGlobalOfCompactFrmt : Error at axis #" << i << " the input big part invalid, end before start !"; + throw INTERP_KERNEL::Exception(oss.str().c_str()); + } + if(partOfBigInAbs[i].first=bigInAbs[i].second) + { + std::ostringstream oss; oss << "MEDCouplingStructuredMesh::ChangeReferenceFromGlobalOfCompactFrmt : Error at axis #" << i << " the part is not included in the big one (start) !"; + throw INTERP_KERNEL::Exception(oss.str().c_str()); + } + partOfBigRelativeToBig[i].first=partOfBigInAbs[i].first-bigInAbs[i].first; + if(partOfBigInAbs[i].secondbigInAbs[i].second) + { + std::ostringstream oss; oss << "MEDCouplingStructuredMesh::ChangeReferenceFromGlobalOfCompactFrmt : Error at axis #" << i << " the part is not included in the big one (end) !"; + throw INTERP_KERNEL::Exception(oss.str().c_str()); + } + partOfBigRelativeToBig[i].second=partOfBigInAbs[i].second-bigInAbs[i].first; + } +} + +/* + * This method is performs the opposite reference modification than explained in ChangeReferenceFromGlobalOfCompactFrmt. + * + * \sa ChangeReferenceFromGlobalOfCompactFrmt + */ +void MEDCouplingStructuredMesh::ChangeReferenceToGlobalOfCompactFrmt(const std::vector< std::pair >& bigInAbs, const std::vector< std::pair >& partOfBigRelativeToBig, std::vector< std::pair >& partOfBigInAbs) +{ + std::size_t dim(bigInAbs.size()); + if(dim!=partOfBigRelativeToBig.size()) + throw INTERP_KERNEL::Exception("MEDCouplingStructuredMesh::ChangeReferenceToGlobalOfCompactFrmt : The size of parts (dimension) must be the same !"); + partOfBigInAbs.resize(dim); + for(std::size_t i=0;ibigInAbs[i].second) + { + std::ostringstream oss; oss << "MEDCouplingStructuredMesh::ChangeReferenceToGlobalOfCompactFrmt : Error at axis #" << i << " the input big part invalid, end before start !"; + throw INTERP_KERNEL::Exception(oss.str().c_str()); + } + if(partOfBigRelativeToBig[i].first<0 || partOfBigRelativeToBig[i].first>=bigInAbs[i].second-bigInAbs[i].first) + { + std::ostringstream oss; oss << "MEDCouplingStructuredMesh::ChangeReferenceToGlobalOfCompactFrmt : Error at axis #" << i << " the start of part is not in the big one !"; + throw INTERP_KERNEL::Exception(oss.str().c_str()); + } + partOfBigInAbs[i].first=partOfBigRelativeToBig[i].first+bigInAbs[i].first; + if(partOfBigRelativeToBig[i].secondbigInAbs[i].second-bigInAbs[i].first) + { + std::ostringstream oss; oss << "MEDCouplingStructuredMesh::ChangeReferenceToGlobalOfCompactFrmt : Error at axis #" << i << " the end of part is not in the big one !"; + throw INTERP_KERNEL::Exception(oss.str().c_str()); + } + partOfBigInAbs[i].second=partOfBigRelativeToBig[i].second+bigInAbs[i].first; + } +} + /*! * This method builds the explicit entity array from the structure in \a st and the range in \a partCompactFormat. * If the range contains invalid values regarding sructure an exception will be thrown. * * \return DataArrayInt * - a new object. - * \sa MEDCouplingStructuredMesh::IsPartStructured, MEDCouplingStructuredMesh::DeduceNumberOfGivenRangeInCompactFrmt, SwitchOnIdsFrom + * \sa MEDCouplingStructuredMesh::IsPartStructured, MEDCouplingStructuredMesh::DeduceNumberOfGivenRangeInCompactFrmt, SwitchOnIdsFrom, ExtractFieldOfBoolFrom */ DataArrayInt *MEDCouplingStructuredMesh::BuildExplicitIdsFrom(const std::vector& st, const std::vector< std::pair >& partCompactFormat) { diff --git a/src/MEDCoupling/MEDCouplingStructuredMesh.hxx b/src/MEDCoupling/MEDCouplingStructuredMesh.hxx index ef24a8546..1175ab37d 100644 --- a/src/MEDCoupling/MEDCouplingStructuredMesh.hxx +++ b/src/MEDCoupling/MEDCouplingStructuredMesh.hxx @@ -74,6 +74,9 @@ namespace ParaMEDMEM MEDCOUPLING_EXPORT static bool IsPartStructured(const int *startIds, const int *stopIds, const std::vector& st, std::vector< std::pair >& partCompactFormat); MEDCOUPLING_EXPORT static std::vector GetDimensionsFromCompactFrmt(const std::vector< std::pair >& partCompactFormat); MEDCOUPLING_EXPORT static void SwitchOnIdsFrom(const std::vector& st, const std::vector< std::pair >& partCompactFormat, std::vector& vectToSwitchOn); + MEDCOUPLING_EXPORT static void ExtractFieldOfBoolFrom(const std::vector& st, const std::vector& fieldOfBool, const std::vector< std::pair >& partCompactFormat, std::vector& fieldOut); + MEDCOUPLING_EXPORT static void ChangeReferenceFromGlobalOfCompactFrmt(const std::vector< std::pair >& bigInAbs, const std::vector< std::pair >& partOfBigInAbs, std::vector< std::pair >& partOfBigRelativeToBig); + MEDCOUPLING_EXPORT static void ChangeReferenceToGlobalOfCompactFrmt(const std::vector< std::pair >& bigInAbs, const std::vector< std::pair >& partOfBigRelativeToBig, std::vector< std::pair >& partOfBigInAbs); MEDCOUPLING_EXPORT static DataArrayInt *BuildExplicitIdsFrom(const std::vector& st, const std::vector< std::pair >& partCompactFormat); MEDCOUPLING_EXPORT static DataArrayInt *Build1GTNodalConnectivity(const int *nodeStBg, const int *nodeStEnd); MEDCOUPLING_EXPORT static DataArrayInt *Build1GTNodalConnectivityOfSubLevelMesh(const int *nodeStBg, const int *nodeStEnd); @@ -81,6 +84,7 @@ namespace ParaMEDMEM 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); + MEDCOUPLING_EXPORT static std::vector< std::vector > ComputeSignaturePerAxisOf(const std::vector& st, const std::vector& crit); private: static int GetNumberOfCellsOfSubLevelMesh(const std::vector& cgs, int mdim); static void GetReverseNodalConnectivity1(const std::vector& ngs, DataArrayInt *revNodal, DataArrayInt *revNodalIndx); @@ -91,9 +95,9 @@ 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 int FindMinimalPartOf1D(const std::vector& st, const std::vector& crit, std::vector< std::pair >& partCompactFormat); + static int FindMinimalPartOf2D(const std::vector& st, const std::vector& crit, std::vector< std::pair >& partCompactFormat); + static int FindMinimalPartOf3D(const std::vector& st, const std::vector& crit, 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]); diff --git a/src/MEDCoupling_Swig/MEDCouplingBasicsTest.py b/src/MEDCoupling_Swig/MEDCouplingBasicsTest.py index c9e409eda..530995096 100644 --- a/src/MEDCoupling_Swig/MEDCouplingBasicsTest.py +++ b/src/MEDCoupling_Swig/MEDCouplingBasicsTest.py @@ -15105,8 +15105,8 @@ class MEDCouplingBasicsTest(unittest.TestCase): fine=DataArrayDouble(3*2*3*4*4*4) ; fine.iota(0) #X=3,Y=2,Z=3 refined by 4 MEDCouplingIMesh.SpreadCoarseToFine(coarse,[5,7,5],fine,[(1,4),(2,4),(1,4)],[4,4,4]) self.assertTrue(fine.isEqual(DataArrayDouble([46.,46.,46.,46.,47.,47.,47.,47.,48.,48.,48.,48.,46.,46.,46.,46.,47.,47.,47.,47.,48.,48.,48.,48.,46.,46.,46.,46.,47.,47.,47.,47.,48.,48.,48.,48.,46.,46.,46.,46.,47.,47.,47.,47.,48.,48.,48.,48.,51.,51.,51.,51.,52.,52.,52.,52.,53.,53.,53.,53.,51.,51.,51.,51.,52.,52.,52.,52.,53.,53.,53.,53.,51.,51.,51.,51.,52.,52.,52.,52.,53.,53.,53.,53.,51.,51.,51.,51.,52.,52.,52.,52.,53.,53.,53.,53.,46.,46.,46.,46.,47.,47.,47.,47.,48.,48.,48.,48.,46.,46.,46.,46.,47.,47.,47.,47.,48.,48.,48.,48.,46.,46.,46.,46.,47.,47.,47.,47.,48.,48.,48.,48.,46.,46.,46.,46.,47.,47.,47.,47.,48.,48.,48.,48.,51.,51.,51.,51.,52.,52.,52.,52.,53.,53.,53.,53.,51.,51.,51.,51.,52.,52.,52.,52.,53.,53.,53.,53.,51.,51.,51.,51.,52.,52.,52.,52.,53.,53.,53.,53.,51.,51.,51.,51.,52.,52.,52.,52.,53.,53.,53.,53.,46.,46.,46.,46.,47.,47.,47.,47.,48.,48.,48.,48.,46.,46.,46.,46.,47.,47.,47.,47.,48.,48.,48.,48.,46.,46.,46.,46.,47.,47.,47.,47.,48.,48.,48.,48.,46.,46.,46.,46.,47.,47.,47.,47.,48.,48.,48.,48.,51.,51.,51.,51.,52.,52.,52.,52.,53.,53.,53.,53.,51.,51.,51.,51.,52.,52.,52.,52.,53.,53.,53.,53.,51.,51.,51.,51.,52.,52.,52.,52.,53.,53.,53.,53.,51.,51.,51.,51.,52.,52.,52.,52.,53.,53.,53.,53.,46.,46.,46.,46.,47.,47.,47.,47.,48.,48.,48.,48.,46.,46.,46.,46.,47.,47.,47.,47.,48.,48.,48.,48.,46.,46.,46.,46.,47.,47.,47.,47.,48.,48.,48.,48.,46.,46.,46.,46.,47.,47.,47.,47.,48.,48.,48.,48.,51.,51.,51.,51.,52.,52.,52.,52.,53.,53.,53.,53.,51.,51.,51.,51.,52.,52.,52.,52.,53.,53.,53.,53.,51.,51.,51.,51.,52.,52.,52.,52.,53.,53.,53.,53.,51.,51.,51.,51.,52.,52.,52.,52.,53.,53.,53.,53.,81.,81.,81.,81.,82.,82.,82.,82.,83.,83.,83.,83.,81.,81.,81.,81.,82.,82.,82.,82.,83.,83.,83.,83.,81.,81.,81.,81.,82.,82.,82.,82.,83.,83.,83.,83.,81.,81.,81.,81.,82.,82.,82.,82.,83.,83.,83.,83.,86.,86.,86.,86.,87.,87.,87.,87.,88.,88.,88.,88.,86.,86.,86.,86.,87.,87.,87.,87.,88.,88.,88.,88.,86.,86.,86.,86.,87.,87.,87.,87.,88.,88.,88.,88.,86.,86.,86.,86.,87.,87.,87.,87.,88.,88.,88.,88.,81.,81.,81.,81.,82.,82.,82.,82.,83.,83.,83.,83.,81.,81.,81.,81.,82.,82.,82.,82.,83.,83.,83.,83.,81.,81.,81.,81.,82.,82.,82.,82.,83.,83.,83.,83.,81.,81.,81.,81.,82.,82.,82.,82.,83.,83.,83.,83.,86.,86.,86.,86.,87.,87.,87.,87.,88.,88.,88.,88.,86.,86.,86.,86.,87.,87.,87.,87.,88.,88.,88.,88.,86.,86.,86.,86.,87.,87.,87.,87.,88.,88.,88.,88.,86.,86.,86.,86.,87.,87.,87.,87.,88.,88.,88.,88.,81.,81.,81.,81.,82.,82.,82.,82.,83.,83.,83.,83.,81.,81.,81.,81.,82.,82.,82.,82.,83.,83.,83.,83.,81.,81.,81.,81.,82.,82.,82.,82.,83.,83.,83.,83.,81.,81.,81.,81.,82.,82.,82.,82.,83.,83.,83.,83.,86.,86.,86.,86.,87.,87.,87.,87.,88.,88.,88.,88.,86.,86.,86.,86.,87.,87.,87.,87.,88.,88.,88.,88.,86.,86.,86.,86.,87.,87.,87.,87.,88.,88.,88.,88.,86.,86.,86.,86.,87.,87.,87.,87.,88.,88.,88.,88.,81.,81.,81.,81.,82.,82.,82.,82.,83.,83.,83.,83.,81.,81.,81.,81.,82.,82.,82.,82.,83.,83.,83.,83.,81.,81.,81.,81.,82.,82.,82.,82.,83.,83.,83.,83.,81.,81.,81.,81.,82.,82.,82.,82.,83.,83.,83.,83.,86.,86.,86.,86.,87.,87.,87.,87.,88.,88.,88.,88.,86.,86.,86.,86.,87.,87.,87.,87.,88.,88.,88.,88.,86.,86.,86.,86.,87.,87.,87.,87.,88.,88.,88.,88.,86.,86.,86.,86.,87.,87.,87.,87.,88.,88.,88.,88.,116.,116.,116.,116.,117.,117.,117.,117.,118.,118.,118.,118.,116.,116.,116.,116.,117.,117.,117.,117.,118.,118.,118.,118.,116.,116.,116.,116.,117.,117.,117.,117.,118.,118.,118.,118.,116.,116.,116.,116.,117.,117.,117.,117.,118.,118.,118.,118.,121.,121.,121.,121.,122.,122.,122.,122.,123.,123.,123.,123.,121.,121.,121.,121.,122.,122.,122.,122.,123.,123.,123.,123.,121.,121.,121.,121.,122.,122.,122.,122.,123.,123.,123.,123.,121.,121.,121.,121.,122.,122.,122.,122.,123.,123.,123.,123.,116.,116.,116.,116.,117.,117.,117.,117.,118.,118.,118.,118.,116.,116.,116.,116.,117.,117.,117.,117.,118.,118.,118.,118.,116.,116.,116.,116.,117.,117.,117.,117.,118.,118.,118.,118.,116.,116.,116.,116.,117.,117.,117.,117.,118.,118.,118.,118.,121.,121.,121.,121.,122.,122.,122.,122.,123.,123.,123.,123.,121.,121.,121.,121.,122.,122.,122.,122.,123.,123.,123.,123.,121.,121.,121.,121.,122.,122.,122.,122.,123.,123.,123.,123.,121.,121.,121.,121.,122.,122.,122.,122.,123.,123.,123.,123.,116.,116.,116.,116.,117.,117.,117.,117.,118.,118.,118.,118.,116.,116.,116.,116.,117.,117.,117.,117.,118.,118.,118.,118.,116.,116.,116.,116.,117.,117.,117.,117.,118.,118.,118.,118.,116.,116.,116.,116.,117.,117.,117.,117.,118.,118.,118.,118.,121.,121.,121.,121.,122.,122.,122.,122.,123.,123.,123.,123.,121.,121.,121.,121.,122.,122.,122.,122.,123.,123.,123.,123.,121.,121.,121.,121.,122.,122.,122.,122.,123.,123.,123.,123.,121.,121.,121.,121.,122.,122.,122.,122.,123.,123.,123.,123.,116.,116.,116.,116.,117.,117.,117.,117.,118.,118.,118.,118.,116.,116.,116.,116.,117.,117.,117.,117.,118.,118.,118.,118.,116.,116.,116.,116.,117.,117.,117.,117.,118.,118.,118.,118.,116.,116.,116.,116.,117.,117.,117.,117.,118.,118.,118.,118.,121.,121.,121.,121.,122.,122.,122.,122.,123.,123.,123.,123.,121.,121.,121.,121.,122.,122.,122.,122.,123.,123.,123.,123.,121.,121.,121.,121.,122.,122.,122.,122.,123.,123.,123.,123.,121.,121.,121.,121.,122.,122.,122.,122.,123.,123.,123.,123.]),1e-12)) - #f=MEDCouplingFieldDouble(ON_CELLS) ; f.setMesh(MEDCouplingIMesh("",3,DataArrayInt([6,8,6]),[0.,0.,0.],DataArrayDouble((1.,1.,1.)))) ; f.setArray(coarse) ; f.setName("tutu") ; f.checkCoherency() ; f.writeVTK("coarse.vti") - #f=MEDCouplingFieldDouble(ON_CELLS) ; f.setMesh(MEDCouplingIMesh("",3,DataArrayInt([13,9,13]),[1.,2.,1.],DataArrayDouble((0.25,0.25,0.25)))) ; f.setArray(fine) ; f.setName("tutu") ; f.checkCoherency() ; f.writeVTK("fine.vti") + f=MEDCouplingFieldDouble(ON_CELLS) ; f.setMesh(MEDCouplingIMesh("",3,DataArrayInt([6,8,6]),[0.,0.,0.],DataArrayDouble((1.,1.,1.)))) ; f.setArray(coarse) ; f.setName("tutu") ; f.checkCoherency() + f=MEDCouplingFieldDouble(ON_CELLS) ; f.setMesh(MEDCouplingIMesh("",3,DataArrayInt([13,9,13]),[1.,2.,1.],DataArrayDouble((0.25,0.25,0.25)))) ; f.setArray(fine) ; f.setName("tutu") ; f.checkCoherency() # 1D coarse=DataArrayDouble(5) ; coarse.iota(0) #X=5 fine=DataArrayDouble(3*4) ; fine.iota(0) #X=3 refined by 4 @@ -15114,6 +15114,30 @@ class MEDCouplingBasicsTest(unittest.TestCase): self.assertTrue(fine.isEqual(DataArrayDouble([1.,1.,1.,1.,2.,2.,2.,2.,3.,3.,3.,3.]),1e-12)) pass + def testAMR4(self): + """This test focuses on MEDCouplingCartesianAMRMesh.createPatchesFromCriterion method. To test it a field containing 0 everywhere except in the annulus (centered on the center of the mesh) value is 1.""" + im=MEDCouplingIMesh("mesh",2,[51,51],[0.,0.],[0.04,0.04]) + b=im.getBarycenterAndOwner() ; b-=[1.,1.] ; b=b.magnitude() + ids=b.getIdsInRange(0.4,0.7) + f=MEDCouplingFieldDouble(ON_CELLS) ; f.setMesh(im) ; f.setName("toto") ; arr=DataArrayDouble(im.getNumberOfCells()) ; arr[:]=0. ; arr[ids]=1. ; f.setArray(arr) + # f.write("test.vti") + amr=MEDCouplingCartesianAMRMesh("mesh",2,[51,51],[0.,0.],[0.04,0.04]) + arr2=DataArrayByte(im.getNumberOfCells()) ; arr2[:]=0 ; arr2[ids]=1 + bso=BoxSplittingOptions() ; bso.setEffeciency(0.8) ; bso.setEffeciencySnd(0.8) ; bso.setMaxCells(1000) ; bso.setMinCellDirection(3) + amr.createPatchesFromCriterion(bso,arr2,[2,2]) + self.assertEqual(18,amr.getNumberOfPatches()) + exp0=[[(8,14),(19,38)],[(19,31),(8,17)],[(19,31),(33,42)],[(10,14),(12,16)],[(9,14),(16,19)],[(14,19),(9,19)],[(14,17),(19,22)],[(14,19),(31,41)],[(36,42),(19,38)],[(14,15),(22,28)],[(14,17),(28,31)],[(31,36),(9,19)],[(33,36),(19,22)],[(31,36),(31,41)],[(36,40),(12,16)],[(36,41),(16,19)],[(35,36),(22,28)],[(33,36),(28,31)]] + for i,bltr in enumerate(exp0): + self.assertEqual(amr[i].getBLTRRange(),bltr) + pass + m=amr.buildMeshFromPatchEnvelop() + self.assertTrue(m.getNodalConnectivity().isEqual(DataArrayInt([1,0,2,3,5,4,6,7,9,8,10,11,13,12,14,15,17,16,18,19,21,20,22,23,25,24,26,27,29,28,30,31,33,32,34,35,37,36,38,39,41,40,42,43,45,44,46,47,49,48,50,51,53,52,54,55,57,56,58,59,61,60,62,63,65,64,66,67,69,68,70,71]))) + self.assertTrue(m.getCoords().isEqualWithoutConsideringStr(DataArrayDouble([0.32,0.76,0.56,0.76,0.32,1.52,0.56,1.52,0.76,0.32,1.24,0.32,0.76,0.68,1.24,0.68,0.76,1.32,1.24,1.32,0.76,1.68,1.24,1.68,0.4,0.48,0.56,0.48,0.4,0.64,0.56,0.64,0.36,0.64,0.56,0.64,0.36,0.76,0.56,0.76,0.56,0.36,0.76,0.36,0.56,0.76,0.76,0.76,0.56,0.76,0.68,0.76,0.56,0.88,0.68,0.88,0.56,1.24,0.76,1.24,0.56,1.64,0.76,1.64,1.44,0.76,1.68,0.76,1.44,1.52,1.68,1.52,0.56,0.88,0.6,0.88,0.56,1.12,0.6,1.12,0.56,1.12,0.68,1.12,0.56,1.24,0.68,1.24,1.24,0.36,1.44,0.36,1.24,0.76,1.44,0.76,1.32,0.76,1.44,0.76,1.32,0.88,1.44,0.88,1.24,1.24,1.44,1.24,1.24,1.64,1.44,1.64,1.44,0.48,1.6,0.48,1.44,0.64,1.6,0.64,1.44,0.64,1.64,0.64,1.44,0.76,1.64,0.76,1.4,0.88,1.44,0.88,1.4,1.12,1.44,1.12,1.32,1.12,1.44,1.12,1.32,1.24,1.44,1.24],72,2),1e-12)) + # + self.assertEqual(MEDCouplingStructuredMesh.ChangeReferenceToGlobalOfCompactFrmt([(8,32),(4,17)],[(0,24),(2,12)]),[(8,32),(6,16)]) + self.assertEqual(MEDCouplingStructuredMesh.ChangeReferenceFromGlobalOfCompactFrmt([(8,32),(4,17)],[(8,32),(6,16)]),[(0,24),(2,12)]) + pass + def setUp(self): pass pass diff --git a/src/MEDCoupling_Swig/MEDCouplingCommon.i b/src/MEDCoupling_Swig/MEDCouplingCommon.i index 4c1ba6391..cc6307736 100644 --- a/src/MEDCoupling_Swig/MEDCouplingCommon.i +++ b/src/MEDCoupling_Swig/MEDCouplingCommon.i @@ -306,6 +306,7 @@ using namespace INTERP_KERNEL; %newobject ParaMEDMEM::MEDCouplingCMesh::clone; %newobject ParaMEDMEM::MEDCouplingCMesh::getCoordsAt; %newobject ParaMEDMEM::MEDCouplingIMesh::New; +%newobject ParaMEDMEM::MEDCouplingIMesh::asSingleCell; %newobject ParaMEDMEM::MEDCouplingIMesh::convertToCartesian; %newobject ParaMEDMEM::MEDCouplingCurveLinearMesh::New; %newobject ParaMEDMEM::MEDCouplingCurveLinearMesh::clone; @@ -317,6 +318,8 @@ using namespace INTERP_KERNEL; %newobject ParaMEDMEM::MEDCouplingCartesianAMRPatch::__getitem__; %newobject ParaMEDMEM::MEDCouplingCartesianAMRMesh::New; %newobject ParaMEDMEM::MEDCouplingCartesianAMRMesh::buildUnstructured; +%newobject ParaMEDMEM::MEDCouplingCartesianAMRMesh::buildMeshFromPatchEnvelop; +%newobject ParaMEDMEM::MEDCouplingCartesianAMRMesh::getImageMesh; %newobject ParaMEDMEM::MEDCouplingCartesianAMRMesh::getGodFather; %newobject ParaMEDMEM::MEDCouplingCartesianAMRMesh::getFather; %newobject ParaMEDMEM::MEDCouplingCartesianAMRMesh::getPatch; @@ -389,6 +392,13 @@ namespace INTERP_KERNEL 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); + %extend + { + std::string __str__() const throw(INTERP_KERNEL::Exception) + { + return self->printOptions(); + } + } }; } @@ -2911,6 +2921,40 @@ namespace ParaMEDMEM PyTuple_SetItem(ret,1,ret1Py); return ret; } + + static PyObject *ChangeReferenceFromGlobalOfCompactFrmt(PyObject *bigInAbs, PyObject *partOfBigInAbs) throw(INTERP_KERNEL::Exception) + { + std::vector< std::pair > param0,param1,ret; + convertPyToVectorPairInt(bigInAbs,param0); + convertPyToVectorPairInt(partOfBigInAbs,param1); + MEDCouplingStructuredMesh::ChangeReferenceFromGlobalOfCompactFrmt(param0,param1,ret); + PyObject *retPy(PyList_New(ret.size())); + for(std::size_t i=0;i > param0,param1,ret; + convertPyToVectorPairInt(bigInAbs,param0); + convertPyToVectorPairInt(partOfBigRelativeToBig,param1); + MEDCouplingStructuredMesh::ChangeReferenceToGlobalOfCompactFrmt(param0,param1,ret); + PyObject *retPy(PyList_New(ret.size())); + for(std::size_t i=0;i& factors) throw(INTERP_KERNEL::Exception); + MEDCouplingIMesh *asSingleCell() const throw(INTERP_KERNEL::Exception); %extend { MEDCouplingIMesh() @@ -4712,8 +4757,10 @@ namespace ParaMEDMEM // int getNumberOfPatches() const throw(INTERP_KERNEL::Exception); MEDCouplingUMesh *buildUnstructured() const throw(INTERP_KERNEL::Exception); + MEDCoupling1SGTUMesh *buildMeshFromPatchEnvelop() const throw(INTERP_KERNEL::Exception); void removePatch(int patchId) throw(INTERP_KERNEL::Exception); void detachFromFather() throw(INTERP_KERNEL::Exception); + void createPatchesFromCriterion(const INTERP_KERNEL::BoxSplittingOptions& bso, const DataArrayByte *criterion, const std::vector& factors) throw(INTERP_KERNEL::Exception); %extend { static MEDCouplingCartesianAMRMesh *New(const std::string& meshName, int spaceDim, PyObject *nodeStrct, PyObject *origin, PyObject *dxyz) throw(INTERP_KERNEL::Exception) @@ -4771,6 +4818,14 @@ namespace ParaMEDMEM return ret; } + MEDCouplingIMesh *getImageMesh() const throw(INTERP_KERNEL::Exception) + { + const MEDCouplingIMesh *ret(self->getImageMesh()); + if(ret) + ret->incrRef(); + return const_cast(ret); + } + MEDCouplingCartesianAMRPatch *__getitem__(int patchId) const throw(INTERP_KERNEL::Exception) { if(patchId==self->getNumberOfPatches()) diff --git a/src/MEDCoupling_Swig/MEDCouplingRemapperTest.py b/src/MEDCoupling_Swig/MEDCouplingRemapperTest.py index 9974af73d..72f92c59e 100644 --- a/src/MEDCoupling_Swig/MEDCouplingRemapperTest.py +++ b/src/MEDCoupling_Swig/MEDCouplingRemapperTest.py @@ -765,11 +765,11 @@ class MEDCouplingBasicsTest(unittest.TestCase): """ This test is the origin of the ref values for MEDCouplingBasicsTest.testAMR2""" coarse=DataArrayDouble(35) ; coarse.iota(0) #X=5,Y=7 fine=DataArrayDouble(3*2*4*4) ; fine.iota(0) #X=3,Y=2 refined by 4 - MEDCouplingIMesh.CondenseFineToCoarse(coarse,[5,7],fine,[(1,4),(2,4)]) + MEDCouplingIMesh.CondenseFineToCoarse(coarse,[5,7],fine,[(1,4),(2,4)],[4,4]) # m=MEDCouplingCartesianAMRMesh("mesh",2,[6,8],[0.,0.],[1.,1.]) trgMesh=m.buildUnstructured() - m.addPatch([(1,4),(2,4)],4) + m.addPatch([(1,4),(2,4)],[4,4]) srcMesh=m[0].getMesh().buildUnstructured() srcField=MEDCouplingFieldDouble(ON_CELLS) fine2=DataArrayDouble(3*2*4*4) ; fine2.iota(0) ; srcField.setArray(fine2)