void MEDCouplingCMesh::setCoordsAt(int i, const DataArrayDouble *arr) throw(INTERP_KERNEL::Exception)
{
+ if(arr)
+ arr->checkNbOfComps(1,"MEDCouplingCMesh::setCoordsAt");
DataArrayDouble **thisArr[3]={&_x_array,&_y_array,&_z_array};
if(i<0 || i>2)
throw INTERP_KERNEL::Exception("Invalid rank specified must be 0 or 1 or 2.");
void MEDCouplingCMesh::setCoords(const DataArrayDouble *coordsX, const DataArrayDouble *coordsY, const DataArrayDouble *coordsZ)
{
+ if(coordsX)
+ coordsX->checkNbOfComps(1,"MEDCouplingCMesh::setCoords : coordsX");
+ if(coordsY)
+ coordsY->checkNbOfComps(1,"MEDCouplingCMesh::setCoords : coordsY");
+ if(coordsZ)
+ coordsZ->checkNbOfComps(1,"MEDCouplingCMesh::setCoords : coordsZ");
if(_x_array)
_x_array->decrRef();
_x_array=const_cast<DataArrayDouble *>(coordsX);
DataArrayInt *MEDCouplingFieldDiscretizationGauss::getOffsetArr(const MEDCouplingMesh *mesh) const
{
int nbOfTuples=mesh->getNumberOfCells();
- DataArrayInt *ret=DataArrayInt::New();
+ MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret=DataArrayInt::New();
ret->alloc(nbOfTuples+1,1);
int *retPtr=ret->getPointer();
const int *start=_discr_per_cell->getConstPointer();
+ if(_discr_per_cell->getNumberOfTuples()!=nbOfTuples)
+ throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getOffsetArr : mismatch between the mesh and the discretization ids array length !");
+ int maxPossible=(int)_loc.size();
retPtr[0]=0;
for(int i=0;i<nbOfTuples;i++,start++)
- retPtr[i+1]=retPtr[i]+_loc[*start].getNumberOfGaussPt();
+ {
+ if(*start>=0 && *start<maxPossible)
+ retPtr[i+1]=retPtr[i]+_loc[*start].getNumberOfGaussPt();
+ else
+ {
+ std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getOffsetArr : At position #" << i << " the locid = " << *start << " whereas it should be in [0," << maxPossible << ") !";
+ throw INTERP_KERNEL::Exception(oss.str().c_str());
+ }
+ }
+ ret->incrRef();
return ret;
}
{
const int *start=_discr_per_cell->getConstPointer();
int nbOfTuples=_discr_per_cell->getNumberOfTuples();
- int *tmp=new int[_loc.size()];
- std::fill(tmp,tmp+_loc.size(),-2);
+ INTERP_KERNEL::AutoPtr<int> tmp=new int[_loc.size()];
+ std::fill((int *)tmp,(int *)tmp+_loc.size(),-2);
for(const int *w=start;w!=start+nbOfTuples;w++)
if(*w>=0)
tmp[*w]=1;
if(tmp[i]!=-2)
tmp[i]=fid++;
if(fid==(int)_loc.size())
- {//no zip needed
- delete [] tmp;
- return;
- }
+ return;
// zip needed
int *start2=_discr_per_cell->getPointer();
for(int *w2=start2;w2!=start2+nbOfTuples;w2++)
- *w2=tmp[*w2];
+ if(*w2>=0)
+ *w2=tmp[*w2];
std::vector<MEDCouplingGaussLocalization> tmpLoc;
for(int i=0;i<(int)_loc.size();i++)
if(tmp[i]!=-2)
tmpLoc.push_back(_loc[tmp[i]]);
- delete [] tmp;
_loc=tmpLoc;
}
{
MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coords=getLocalizationOfDiscValues(mesh);
int nbOfPts=coords->getNumberOfTuples();
- int dimension=coords->getNumberOfComponents();
+ //int dimension=coords->getNumberOfComponents();
MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrix=coords->buildEuclidianDistanceDenseMatrix();
operateOnDenseMatrix(mesh->getSpaceDimension(),nbOfPts*nbOfPts,matrix->getPointer());
// Drift
}
template<int SPACEDIM>
-void DataArrayDouble::findTupleIdsNearTuplesAlg(const BBTree<SPACEDIM,int>& myTree, const double *pos, int nbOfTuples, double eps,
- std::vector<int>& c, std::vector<int>& cI) const
+void DataArrayDouble::FindTupleIdsNearTuplesAlg(const BBTree<SPACEDIM,int>& myTree, const double *pos, int nbOfTuples, double eps,
+ std::vector<int>& c, std::vector<int>& cI)
{
- const double *coordsPtr=getConstPointer();
for(int i=0;i<nbOfTuples;i++)
{
std::vector<int> intersectingElems;
case 3:
{
BBTree<3,int> myTree(bbox->getConstPointer(),0,0,getNumberOfTuples(),eps/10);
- findTupleIdsNearTuplesAlg<3>(myTree,other->getConstPointer(),nbOfTuplesOther,eps,c,cI);
+ FindTupleIdsNearTuplesAlg<3>(myTree,other->getConstPointer(),nbOfTuplesOther,eps,c,cI);
break;
}
case 2:
{
BBTree<2,int> myTree(bbox->getConstPointer(),0,0,getNumberOfTuples(),eps/10);
- findTupleIdsNearTuplesAlg<2>(myTree,other->getConstPointer(),nbOfTuplesOther,eps,c,cI);
+ FindTupleIdsNearTuplesAlg<2>(myTree,other->getConstPointer(),nbOfTuplesOther,eps,c,cI);
break;
}
case 1:
{
BBTree<1,int> myTree(bbox->getConstPointer(),0,0,getNumberOfTuples(),eps/10);
- findTupleIdsNearTuplesAlg<1>(myTree,other->getConstPointer(),nbOfTuplesOther,eps,c,cI);
+ FindTupleIdsNearTuplesAlg<1>(myTree,other->getConstPointer(),nbOfTuplesOther,eps,c,cI);
break;
}
default:
MEDCOUPLING_EXPORT int getNumberOfTuples() const { return _nb_of_tuples; }
MEDCOUPLING_EXPORT int getNbOfElems() const { return ((int)_info_on_compo.size())*_nb_of_tuples; }
MEDCOUPLING_EXPORT void checkNbOfTuples(int nbOfTuples, const char *msg) const throw(INTERP_KERNEL::Exception);
- MEDCOUPLING_EXPORT
-void checkNbOfComps(int nbOfCompo, const char *msg) const throw(INTERP_KERNEL::Exception);
+ MEDCOUPLING_EXPORT void checkNbOfComps(int nbOfCompo, const char *msg) const throw(INTERP_KERNEL::Exception);
MEDCOUPLING_EXPORT void checkNbOfTuplesAndComp(const DataArray& other, const char *msg) const throw(INTERP_KERNEL::Exception);
MEDCOUPLING_EXPORT void checkNbOfTuplesAndComp(int nbOfTuples, int nbOfCompo, const char *msg) const throw(INTERP_KERNEL::Exception);
MEDCOUPLING_EXPORT void checkNbOfElems(int nbOfElems, const char *msg) const throw(INTERP_KERNEL::Exception);
template<int SPACEDIM>
void findCommonTuplesAlg(const double *bbox, int nbNodes, int limitNodeId, double prec, std::vector<int>& c, std::vector<int>& cI) const;
template<int SPACEDIM>
- void findTupleIdsNearTuplesAlg(const BBTree<SPACEDIM,int>& myTree, const double *pos, int nbOfTuples, double eps,
- std::vector<int>& c, std::vector<int>& cI) const;
+ static void FindTupleIdsNearTuplesAlg(const BBTree<SPACEDIM,int>& myTree, const double *pos, int nbOfTuples, double eps,
+ std::vector<int>& c, std::vector<int>& cI);
private:
DataArrayDouble() { }
private:
throw INTERP_KERNEL::Exception("MEDCouplingUMesh::simplifyPolyhedra : works on meshdimension 3 and spaceDimension 3 !");
MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coords=getCoords()->deepCpy();
coords->recenterForMaxPrecision(eps);
- const double *coordsPtr=coords->getConstPointer();
//
int nbOfCells=getNumberOfCells();
const int *conn=_nodal_connec->getConstPointer();
}
int thisNbCells=getNumberOfCells();
int otherNbCells=other->getNumberOfCells();
- int nbOfCells=mesh->getNumberOfCells();
MEDCouplingAutoRefCountObjectPtr<DataArrayInt> arr2=DataArrayInt::New();
arr2->alloc(otherNbCells,1);
arr2->fillWithZero();
}
int nbOfCells=getNumberOfCells();
bool easyAssign=true;
- const int *conn=_nodal_connec->getConstPointer();
const int *connI=_nodal_connec_index->getConstPointer();
- const int *connOther=otherOnSameCoordsThanThis._nodal_connec->getConstPointer();
const int *connIOther=otherOnSameCoordsThanThis._nodal_connec_index->getConstPointer();
for(const int *it=cellIdsBg;it!=cellIdsEnd && easyAssign;it++,connIOther++)
{
}
int nbOfCells=getNumberOfCells();
bool easyAssign=true;
- const int *conn=_nodal_connec->getConstPointer();
const int *connI=_nodal_connec_index->getConstPointer();
- const int *connOther=otherOnSameCoordsThanThis._nodal_connec->getConstPointer();
const int *connIOther=otherOnSameCoordsThanThis._nodal_connec_index->getConstPointer();
int it=start;
for(int i=0;i<nbOfCellsToModify && easyAssign;i++,it+=step,connIOther++)
{
const int *ptI=_nodal_connec_index->getConstPointer();
const int *pt=_nodal_connec->getConstPointer();
- return (INTERP_KERNEL::NormalizedCellType) pt[ptI[cellId]];
+ if(cellId>=0 && cellId<_nodal_connec_index->getNbOfElems()-1)
+ return (INTERP_KERNEL::NormalizedCellType) pt[ptI[cellId]];
+ else
+ {
+ std::ostringstream oss; oss << "MEDCouplingUMesh::getTypeOfCell : Requesting type of cell #" << cellId << " but it should be in [0," << _nodal_connec_index->getNbOfElems()-1 << ") !";
+ throw INTERP_KERNEL::Exception(oss.str().c_str());
+ }
}
/*!
}
/*!
- * This method makes the hypothesis that 'this' is sorted by type. If not an exception will be thrown.
- * This method is the opposite of MEDCouplingUMesh::checkTypeConsistencyAndContig method. Given a list of cells in 'profile' it returns a list of profiles sorted by geo type.
- * This method has 1 input 'profile' and 2 outputs 'code' and 'idsPerType'.
- * @throw if 'profile' has not exactly one component. It throws too, if 'profile' contains some values not in [0,getNumberOfCells()) or if 'this' is not fully defined
+ * This method makes the hypothesis that \at this is sorted by type. If not an exception will be thrown.
+ * This method is the opposite of MEDCouplingUMesh::checkTypeConsistencyAndContig method. Given a list of cells in \a profile it returns a list of sub-profiles sorted by geo type.
+ * The result is put in the array \a idsPerType. In the returned parameter \a code, foreach i \a code[3*i+2] refers (if different from -1) to a location into the \a idsPerType.
+ * This method has 1 input \a profile and 3 outputs \a code' \a idsInPflPerType and \a idsPerType.
+ *
+ * @param [out] code is a vector of size 3*n where n is the number of different geometric type in \a this \b reduced to the profile \a profile. \a code has exactly the same semantic than in MEDCouplingUMesh::checkTypeConsistencyAndContig method.
+ * @param [out] idsInPflPerType is a vector of size of different geometric type in the subpart defined by \a profile of \a this ( equal to \a code.size()/3). For each i,
+ * \a idsInPflPerType[i] stores the tuple ids in \a profile that correspond to the geometric type code[3*i+0]
+ * @param [out] idsPerType is a vector of size of different sub profiles needed to be defined to represent the profile \a profile for a given geometric type.
+ * This vector can be empty in case of all geometric type cells are fully covered in ascending in the given input \a profile.
+ * @throw if \a profile has not exactly one component. It throws too, if \a profile contains some values not in [0,getNumberOfCells()) or if 'this' is not fully defined
*/
void MEDCouplingUMesh::splitProfilePerType(const DataArrayInt *profile, std::vector<int>& code, std::vector<DataArrayInt *>& idsInPflPerType, std::vector<DataArrayInt *>& idsPerType) const throw(INTERP_KERNEL::Exception)
{
/*!
* This method reduced number of cells of this by keeping cells whose type is different from 'type' and if type=='type'
+ * This method \b works \b for mesh sorted by type.
* cells whose ids is in 'idsPerGeoType' array.
* This method conserves coords and name of mesh.
*/
MEDCouplingUMesh *MEDCouplingUMesh::keepSpecifiedCells(INTERP_KERNEL::NormalizedCellType type, const int *idsPerGeoTypeBg, const int *idsPerGeoTypeEnd) const
{
- std::vector<int> idsTokeep;
- int nbOfCells=getNumberOfCells();
- int j=0;
- for(int i=0;i<nbOfCells;i++)
- if(getTypeOfCell(i)!=type)
- idsTokeep.push_back(i);
- else
+ std::vector<int> code=getDistributionOfTypes();
+ std::size_t nOfTypesInThis=code.size()/3;
+ int sz=0,szOfType=0;
+ for(std::size_t i=0;i<nOfTypesInThis;i++)
+ {
+ if(code[3*i]!=type)
+ sz+=code[3*i+1];
+ else
+ szOfType=code[3*i+1];
+ }
+ for(const int *work=idsPerGeoTypeBg;work!=idsPerGeoTypeEnd;work++)
+ if(*work<0 || *work>=szOfType)
{
- if(std::find(idsPerGeoTypeBg,idsPerGeoTypeEnd,j)!=idsPerGeoTypeEnd)
- idsTokeep.push_back(i);
- j++;
+ std::ostringstream oss; oss << "MEDCouplingUMesh::keepSpecifiedCells : Request on type " << type << " at place #" << std::distance(idsPerGeoTypeBg,work) << " value " << *work;
+ oss << ". It should be in [0," << szOfType << ") !";
+ throw INTERP_KERNEL::Exception(oss.str().c_str());
}
- MEDCouplingPointSet *ret=buildPartOfMySelf(&idsTokeep[0],&idsTokeep[0]+idsTokeep.size(),true);
- MEDCouplingUMesh *ret2=dynamic_cast<MEDCouplingUMesh *>(ret);
- if(!ret2)
+ MEDCouplingAutoRefCountObjectPtr<DataArrayInt> idsTokeep=DataArrayInt::New(); idsTokeep->alloc(sz+(int)std::distance(idsPerGeoTypeBg,idsPerGeoTypeEnd),1);
+ int *idsPtr=idsTokeep->getPointer();
+ int offset=0;
+ for(std::size_t i=0;i<nOfTypesInThis;i++)
{
- ret->decrRef();
- return 0;
+ if(code[3*i]!=type)
+ for(int j=0;j<code[3*i+1];j++)
+ *idsPtr++=offset+j;
+ else
+ idsPtr=std::transform(idsPerGeoTypeBg,idsPerGeoTypeEnd,idsPtr,std::bind2nd(std::plus<int>(),offset));
+ offset+=code[3*i+1];
}
- ret2->copyTinyInfoFrom(this);
- return ret2;
+ MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> ret=static_cast<MEDCouplingUMesh *>(buildPartOfMySelf(idsTokeep->begin(),idsTokeep->end(),true));
+ ret->copyTinyInfoFrom(this);
+ ret->incrRef();
+ return ret;
}
/*!
MEDCouplingAutoRefCountObjectPtr<DataArrayInt> connectivity=DataArrayInt::New(); connectivity->alloc(_nodal_connec->getNumberOfTuples()-nbOfCells,1);
int *w1=faceoffsets->getPointer(),*w2=types->getPointer(),*w3=offsets->getPointer(),*w4=connectivity->getPointer();
int szFaceOffsets=0,szConn=0;
- for(int i=0;i<nbOfCells;i++,w1++,w2++,*w3++)
+ for(int i=0;i<nbOfCells;i++,w1++,w2++,w3++)
{
*w2=cPtr[cIPtr[i]];
if((INTERP_KERNEL::NormalizedCellType)cPtr[cIPtr[i]]!=INTERP_KERNEL::NORM_POLYHED)
{
int nbFaces=std::count(cPtr+cIPtr[i]+1,cPtr+cIPtr[i+1],-1)+1;
*w1++=nbFaces;
- const int *w4=cPtr+cIPtr[i]+1,*w5=0;
+ const int *w6=cPtr+cIPtr[i]+1,*w5=0;
for(int j=0;j<nbFaces;j++)
{
- w5=std::find(w4,cPtr+cIPtr[i+1],-1);
- *w1++=(int)std::distance(w4,w5);
- w1=std::copy(w4,w5,w1);
- w4=w5+1;
+ w5=std::find(w6,cPtr+cIPtr[i+1],-1);
+ *w1++=(int)std::distance(w6,w5);
+ w1=std::copy(w6,w5,w1);
+ w6=w5+1;
}
}
faces->writeVTK(ofs,8,"Int32","faces");
return true;
}
}
+ else
+ throw INTERP_KERNEL::Exception("MEDCouplingUMesh::BuildConvexEnvelopOf2DCellJarvis : invalid 2D cell connectivity !");
}
else
throw INTERP_KERNEL::Exception("MEDCouplingUMesh::BuildConvexEnvelopOf2DCellJarvis : invalid 2D cell connectivity !");
int spaceDim=getSpaceDimension();
if(mdim!=spaceDim)
throw INTERP_KERNEL::Exception("MEDCouplingUMesh::buildSpreadZonesWithPoly : meshdimension and spacedimension do not match !");
- int nbCells=getNumberOfCells();
std::vector<DataArrayInt *> partition=partitionBySpreadZone();
std::vector< MEDCouplingAutoRefCountObjectPtr<DataArrayInt> > partitionAuto; partitionAuto.reserve(partition.size());
std::copy(partition.begin(),partition.end(),std::back_insert_iterator<std::vector< MEDCouplingAutoRefCountObjectPtr<DataArrayInt> > >(partitionAuto));
DataArrayDouble *arr4=DataArrayDouble::New();
arr4->alloc(7,2);
double *tmp2=arr4->getPointer();
- const int arr4Ref[14]={0.8,10.8,1.9,11.9,2.1,12.1,3.2,13.2,4.3,14.3,5.4,15.4,6.5,16.5};
+ const double arr4Ref[14]={0.8,10.8,1.9,11.9,2.1,12.1,3.2,13.2,4.3,14.3,5.4,15.4,6.5,16.5};
std::copy(arr4Ref,arr4Ref+14,tmp2);
CPPUNIT_ASSERT_EQUAL(7,arr4->getNumberOfTuples());
CPPUNIT_ASSERT_EQUAL(2,arr4->getNumberOfComponents());
CPPUNIT_ASSERT_EQUAL(5,f->getNumberOfMeshPlacesExpected());
CPPUNIT_ASSERT_EQUAL(0,f->getNbOfGaussLocalization());
f->setGaussLocalizationOnType(INTERP_KERNEL::NORM_TRI3,_refCoo1,_gsCoo1,_wg1);
+ f->setGaussLocalizationOnType(INTERP_KERNEL::NORM_TRI3,_refCoo1,_gsCoo1,_wg1); // not a bug only to check that it works well
CPPUNIT_ASSERT_THROW(f->setGaussLocalizationOnType(INTERP_KERNEL::NORM_QUAD4,_refCoo1,_gsCoo1,_wg1),INTERP_KERNEL::Exception);
CPPUNIT_ASSERT_EQUAL(1,f->getNbOfGaussLocalization());
const double refCoo2[8]={ 0.,0., 1.,0., 1.,1., 0.,1. };
DataArrayDouble *d=DataArrayDouble::New();
d->alloc(13,1);
d->iota();
- MEDCouplingCMesh *e=MEDCouplingCMesh::New();
- e->setCoordsAt(0,d);
- MEDCouplingUMesh *f=e->buildUnstructured();
+ MEDCouplingCMesh *ee=MEDCouplingCMesh::New();
+ ee->setCoordsAt(0,d);
+ MEDCouplingUMesh *f=ee->buildUnstructured();
DataArrayDouble *g=f->getCoords()->applyFunc(2,"3.5*IVec+x/6*3.14159265359*JVec");
CPPUNIT_ASSERT_THROW(f->getCoords()->applyFunc(2,"3.5*IVec+x/6*3.14159265359*KVec"),INTERP_KERNEL::Exception); // KVec refers to component #2 and there is only 2 components !
DataArrayDouble *h=g->fromPolarToCart();
h->decrRef();
g->decrRef();
f->decrRef();
- e->decrRef();
+ ee->decrRef();
d->decrRef();
c->decrRef();
b->decrRef();
self.assertEqual(5,f.getNumberOfMeshPlacesExpected());
self.assertEqual(0,f.getNbOfGaussLocalization());
f.setGaussLocalizationOnType(NORM_TRI3,_refCoo1,_gsCoo1,_wg1);
+ f.setGaussLocalizationOnType(NORM_TRI3,_refCoo1,_gsCoo1,_wg1); # not a bug only to check that it works well
self.assertRaises(InterpKernelException,f.setGaussLocalizationOnType,NORM_QUAD4,_refCoo1,_gsCoo1,_wg1)
self.assertEqual(1,f.getNbOfGaussLocalization());
refCoo2=[ 0.,0., 1.,0., 1.,1., 0.,1. ]
fc.getArray()
pass
+ def testSwigRotate(self):
+ d=DataArrayDouble([1.,2.,3.,4.,6.,5.],2,3)
+ MEDCouplingPointSet.Rotate3DAlg([0.,0.,0.],[0.,1.,0.],1.5707963267948966,d)
+ self.assertTrue(d.isEqual(DataArrayDouble([3.,2.,-1.,5.,6.,-4.],2,3),1e-12))
+ d=DataArrayDouble([1.,2.,3.,4.,6.,5.],3,2)
+ MEDCouplingPointSet.Rotate2DAlg([0.,0.],1.5707963267948966,d)
+ self.assertTrue(d.isEqual(DataArrayDouble([-2.,1.,-4.,3.,-5.,6.],3,2),1e-12))
+ pass
+
+ def testSwigCMeshProtection(self):
+ cm=MEDCouplingCMesh()
+ self.assertRaises(InterpKernelException,cm.setCoordsAt,0,DataArrayDouble([4.,4.5,6.,7.],2,2))
+ self.assertRaises(InterpKernelException,cm.setCoords,DataArrayDouble([4.,4.5,6.,7.],2,2))
+ pass
+
def setUp(self):
pass
pass
for(int i=0;i<sz;i++)
PyList_SetItem(coords,i,PyFloat_FromDouble(coo[i]));
}
+
+ static void Rotate2DAlg(PyObject *center, double angle, PyObject *coords) throw(INTERP_KERNEL::Exception)
+ {
+ int sz,sz2;
+ INTERP_KERNEL::AutoPtr<double> c=convertPyToNewDblArr2(center,&sz);
+ int sw,nbNodes=0;
+ double val0; ParaMEDMEM::DataArrayDouble *val1=0; ParaMEDMEM::DataArrayDoubleTuple *val2=0;
+ std::vector<double> val3;
+ const double *coo=convertObjToPossibleCpp5_Safe2(coords,sw,val0,val1,val2,val3,
+ "Rotate2DAlg",2,true,nbNodes);
+ if(sw!=2 && sw!=3)
+ throw INTERP_KERNEL::Exception("Invalid call to MEDCouplingPointSet::Rotate2DAlg : try another overload method !");
+ ParaMEDMEM::MEDCouplingPointSet::Rotate2DAlg(c,angle,nbNodes,const_cast<double *>(coo));
+ }
+
static void Rotate3DAlg(PyObject *center, PyObject *vect, double angle, int nbNodes, PyObject *coords) throw(INTERP_KERNEL::Exception)
{
int sz,sz2;
INTERP_KERNEL::AutoPtr<double> c=convertPyToNewDblArr2(center,&sz);
INTERP_KERNEL::AutoPtr<double> coo=convertPyToNewDblArr2(coords,&sz);
- double *v=convertPyToNewDblArr2(vect,&sz2);
+ INTERP_KERNEL::AutoPtr<double> v=convertPyToNewDblArr2(vect,&sz2);
ParaMEDMEM::MEDCouplingPointSet::Rotate3DAlg(c,v,angle,nbNodes,coo);
for(int i=0;i<sz;i++)
PyList_SetItem(coords,i,PyFloat_FromDouble(coo[i]));
}
+
+ static void Rotate3DAlg(PyObject *center, PyObject *vect, double angle, PyObject *coords) throw(INTERP_KERNEL::Exception)
+ {
+ int sz,sz2;
+ INTERP_KERNEL::AutoPtr<double> c=convertPyToNewDblArr2(center,&sz);
+ int sw,nbNodes=0;
+ double val0; ParaMEDMEM::DataArrayDouble *val1=0; ParaMEDMEM::DataArrayDoubleTuple *val2=0;
+ std::vector<double> val3;
+ const double *coo=convertObjToPossibleCpp5_Safe2(coords,sw,val0,val1,val2,val3,
+ "Rotate3DAlg",3,true,nbNodes);
+ if(sw!=2 && sw!=3)
+ throw INTERP_KERNEL::Exception("Invalid call to MEDCouplingPointSet::Rotate3DAlg : try another overload method !");
+ INTERP_KERNEL::AutoPtr<double> v=convertPyToNewDblArr2(vect,&sz2);
+ ParaMEDMEM::MEDCouplingPointSet::Rotate3DAlg(c,v,angle,nbNodes,const_cast<double *>(coo));
+ }
}
};
for(int i=0;i<size;i++)
{
PyObject *obj=PyTuple_GetItem(pyLi,i);
- void *argp;
int status=SWIG_ConvertPtr(obj,&argp,ty,0|0);
if(!SWIG_IsOK(status))
{