1 // Copyright (C) 2007-2012 CEA/DEN, EDF R&D
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License.
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
11 // Lesser General Public License for more details.
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
19 // Author : Anthony Geay (CEA/DEN)
21 #include "MEDCouplingFieldDiscretization.hxx"
22 #include "MEDCouplingCMesh.hxx"
23 #include "MEDCouplingUMesh.hxx"
24 #include "MEDCouplingFieldDouble.hxx"
25 #include "MEDCouplingAutoRefCountObjectPtr.hxx"
27 #include "CellModel.hxx"
28 #include "InterpolationUtils.hxx"
29 #include "InterpKernelAutoPtr.hxx"
30 #include "InterpKernelGaussCoords.hxx"
31 #include "InterpKernelMatrixTools.hxx"
40 using namespace ParaMEDMEM;
42 const double MEDCouplingFieldDiscretization::DFLT_PRECISION=1.e-12;
44 const char MEDCouplingFieldDiscretizationP0::REPR[]="P0";
46 const TypeOfField MEDCouplingFieldDiscretizationP0::TYPE=ON_CELLS;
48 const char MEDCouplingFieldDiscretizationP1::REPR[]="P1";
50 const TypeOfField MEDCouplingFieldDiscretizationP1::TYPE=ON_NODES;
52 const int MEDCouplingFieldDiscretizationPerCell::DFT_INVALID_LOCID_VALUE=-1;
54 const char MEDCouplingFieldDiscretizationGauss::REPR[]="GAUSS";
56 const TypeOfField MEDCouplingFieldDiscretizationGauss::TYPE=ON_GAUSS_PT;
58 const char MEDCouplingFieldDiscretizationGaussNE::REPR[]="GSSNE";
60 const TypeOfField MEDCouplingFieldDiscretizationGaussNE::TYPE=ON_GAUSS_NE;
62 const char MEDCouplingFieldDiscretizationKriging::REPR[]="KRIGING";
64 const TypeOfField MEDCouplingFieldDiscretizationKriging::TYPE=ON_NODES_KR;
66 MEDCouplingFieldDiscretization::MEDCouplingFieldDiscretization():_precision(DFLT_PRECISION)
70 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretization::New(TypeOfField type)
74 case MEDCouplingFieldDiscretizationP0::TYPE:
75 return new MEDCouplingFieldDiscretizationP0;
76 case MEDCouplingFieldDiscretizationP1::TYPE:
77 return new MEDCouplingFieldDiscretizationP1;
78 case MEDCouplingFieldDiscretizationGauss::TYPE:
79 return new MEDCouplingFieldDiscretizationGauss;
80 case MEDCouplingFieldDiscretizationGaussNE::TYPE:
81 return new MEDCouplingFieldDiscretizationGaussNE;
82 case MEDCouplingFieldDiscretizationKriging::TYPE:
83 return new MEDCouplingFieldDiscretizationKriging;
85 throw INTERP_KERNEL::Exception("Choosen discretization is not implemented yet.");
89 TypeOfField MEDCouplingFieldDiscretization::getTypeOfFieldFromStringRepr(const char *repr) throw(INTERP_KERNEL::Exception)
91 std::string reprCpp(repr);
92 if(reprCpp==MEDCouplingFieldDiscretizationP0::REPR)
93 return MEDCouplingFieldDiscretizationP0::TYPE;
94 if(reprCpp==MEDCouplingFieldDiscretizationP1::REPR)
95 return MEDCouplingFieldDiscretizationP1::TYPE;
96 if(reprCpp==MEDCouplingFieldDiscretizationGauss::REPR)
97 return MEDCouplingFieldDiscretizationGauss::TYPE;
98 if(reprCpp==MEDCouplingFieldDiscretizationGaussNE::REPR)
99 return MEDCouplingFieldDiscretizationGaussNE::TYPE;
100 if(reprCpp==MEDCouplingFieldDiscretizationKriging::REPR)
101 return MEDCouplingFieldDiscretizationKriging::TYPE;
102 throw INTERP_KERNEL::Exception("Representation does not match with any field discretization !");
105 bool MEDCouplingFieldDiscretization::isEqual(const MEDCouplingFieldDiscretization *other, double eps) const
108 return isEqualIfNotWhy(other,eps,reason);
111 bool MEDCouplingFieldDiscretization::isEqualWithoutConsideringStr(const MEDCouplingFieldDiscretization *other, double eps) const
113 return isEqual(other,eps);
117 * For all field discretization excepted GaussPts the [ \a startCellIds, \a endCellIds ) has no impact on the cloned instance.
119 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretization::clonePart(const int *startCellIds, const int *endCellIds) const
125 * Excepted for MEDCouplingFieldDiscretizationPerCell no underlying TimeLabel object : nothing to do in generally.
127 void MEDCouplingFieldDiscretization::updateTime() const
132 * Computes normL1 of DataArrayDouble instance arr.
133 * @param res output parameter expected to be of size arr->getNumberOfComponents();
134 * @throw when the field discretization fails on getMeasure fields (gauss points for example)
136 void MEDCouplingFieldDiscretization::normL1(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, double *res) const throw(INTERP_KERNEL::Exception)
138 MEDCouplingFieldDouble *vol=getMeasureField(mesh,true);
139 int nbOfCompo=arr->getNumberOfComponents();
140 int nbOfElems=getNumberOfTuples(mesh);
141 std::fill(res,res+nbOfCompo,0.);
142 const double *arrPtr=arr->getConstPointer();
143 const double *volPtr=vol->getArray()->getConstPointer();
145 for(int i=0;i<nbOfElems;i++)
147 double v=fabs(volPtr[i]);
148 for(int j=0;j<nbOfCompo;j++)
149 res[j]+=fabs(arrPtr[i*nbOfCompo+j])*v;
152 std::transform(res,res+nbOfCompo,res,std::bind2nd(std::multiplies<double>(),1./deno));
157 * Computes normL2 of DataArrayDouble instance arr.
158 * @param res output parameter expected to be of size arr->getNumberOfComponents();
159 * @throw when the field discretization fails on getMeasure fields (gauss points for example)
161 void MEDCouplingFieldDiscretization::normL2(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, double *res) const throw(INTERP_KERNEL::Exception)
163 MEDCouplingFieldDouble *vol=getMeasureField(mesh,true);
164 int nbOfCompo=arr->getNumberOfComponents();
165 int nbOfElems=getNumberOfTuples(mesh);
166 std::fill(res,res+nbOfCompo,0.);
167 const double *arrPtr=arr->getConstPointer();
168 const double *volPtr=vol->getArray()->getConstPointer();
170 for(int i=0;i<nbOfElems;i++)
172 double v=fabs(volPtr[i]);
173 for(int j=0;j<nbOfCompo;j++)
174 res[j]+=arrPtr[i*nbOfCompo+j]*arrPtr[i*nbOfCompo+j]*v;
177 std::transform(res,res+nbOfCompo,res,std::bind2nd(std::multiplies<double>(),1./deno));
178 std::transform(res,res+nbOfCompo,res,std::ptr_fun<double,double>(std::sqrt));
183 * Computes integral of DataArrayDouble instance arr.
184 * @param res output parameter expected to be of size arr->getNumberOfComponents();
185 * @throw when the field discretization fails on getMeasure fields (gauss points for example)
187 void MEDCouplingFieldDiscretization::integral(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, bool isWAbs, double *res) const throw(INTERP_KERNEL::Exception)
189 MEDCouplingFieldDouble *vol=getMeasureField(mesh,isWAbs);
190 int nbOfCompo=arr->getNumberOfComponents();
191 int nbOfElems=getNumberOfTuples(mesh);
192 std::fill(res,res+nbOfCompo,0.);
193 const double *arrPtr=arr->getConstPointer();
194 const double *volPtr=vol->getArray()->getConstPointer();
195 double *tmp=new double[nbOfCompo];
196 for (int i=0;i<nbOfElems;i++)
198 std::transform(arrPtr+i*nbOfCompo,arrPtr+(i+1)*nbOfCompo,tmp,std::bind2nd(std::multiplies<double>(),volPtr[i]));
199 std::transform(tmp,tmp+nbOfCompo,res,res,std::plus<double>());
205 void MEDCouplingFieldDiscretization::getSerializationIntArray(DataArrayInt *& arr) const
213 void MEDCouplingFieldDiscretization::getTinySerializationIntInformation(std::vector<int>& tinyInfo) const
220 void MEDCouplingFieldDiscretization::getTinySerializationDbleInformation(std::vector<double>& tinyInfo) const
224 void MEDCouplingFieldDiscretization::resizeForUnserialization(const std::vector<int>& tinyInfo, DataArrayInt *& arr)
232 void MEDCouplingFieldDiscretization::finishUnserialization(const std::vector<double>& tinyInfo)
237 * This method is typically the first step of renumbering. The implementation is empty it is not a bug only gauss is impacted
238 * virtualy by this method.
240 void MEDCouplingFieldDiscretization::renumberCells(const int *old2NewBg, bool check) throw(INTERP_KERNEL::Exception)
244 double MEDCouplingFieldDiscretization::getIJK(const MEDCouplingMesh *mesh, const DataArrayDouble *da,
245 int cellId, int nodeIdInCell, int compoId) const throw(INTERP_KERNEL::Exception)
247 throw INTERP_KERNEL::Exception("getIJK Invalid ! only for GaussPoint and GaussNE discretizations !");
250 void MEDCouplingFieldDiscretization::setGaussLocalizationOnType(const MEDCouplingMesh *m, INTERP_KERNEL::NormalizedCellType type, const std::vector<double>& refCoo,
251 const std::vector<double>& gsCoo, const std::vector<double>& wg) throw(INTERP_KERNEL::Exception)
253 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
256 void MEDCouplingFieldDiscretization::setGaussLocalizationOnCells(const MEDCouplingMesh *m, const int *begin, const int *end, const std::vector<double>& refCoo,
257 const std::vector<double>& gsCoo, const std::vector<double>& wg) throw(INTERP_KERNEL::Exception)
259 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
262 void MEDCouplingFieldDiscretization::clearGaussLocalizations() throw(INTERP_KERNEL::Exception)
264 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
267 MEDCouplingGaussLocalization& MEDCouplingFieldDiscretization::getGaussLocalization(int locId) throw(INTERP_KERNEL::Exception)
269 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
272 const MEDCouplingGaussLocalization& MEDCouplingFieldDiscretization::getGaussLocalization(int locId) const throw(INTERP_KERNEL::Exception)
274 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
277 int MEDCouplingFieldDiscretization::getNbOfGaussLocalization() const throw(INTERP_KERNEL::Exception)
279 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
282 int MEDCouplingFieldDiscretization::getGaussLocalizationIdOfOneCell(int cellId) const throw(INTERP_KERNEL::Exception)
284 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
287 int MEDCouplingFieldDiscretization::getGaussLocalizationIdOfOneType(INTERP_KERNEL::NormalizedCellType type) const throw(INTERP_KERNEL::Exception)
289 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
292 std::set<int> MEDCouplingFieldDiscretization::getGaussLocalizationIdsOfOneType(INTERP_KERNEL::NormalizedCellType type) const throw(INTERP_KERNEL::Exception)
294 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
297 void MEDCouplingFieldDiscretization::getCellIdsHavingGaussLocalization(int locId, std::vector<int>& cellIds) const throw(INTERP_KERNEL::Exception)
299 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
302 void MEDCouplingFieldDiscretization::RenumberEntitiesFromO2NArr(double eps, const int *old2NewPtr, int newNbOfEntity, DataArrayDouble *arr, const char *msg)
304 int oldNbOfElems=arr->getNumberOfTuples();
305 int nbOfComp=arr->getNumberOfComponents();
306 int newNbOfTuples=newNbOfEntity;
307 DataArrayDouble *arrCpy=arr->deepCpy();
308 const double *ptSrc=arrCpy->getConstPointer();
309 arr->reAlloc(newNbOfTuples);
310 double *ptToFill=arr->getPointer();
311 std::fill(ptToFill,ptToFill+nbOfComp*newNbOfTuples,std::numeric_limits<double>::max());
312 INTERP_KERNEL::AutoPtr<double> tmp=new double[nbOfComp];
313 for(int i=0;i<oldNbOfElems;i++)
315 int newNb=old2NewPtr[i];
316 if(newNb>=0)//if newNb<0 the node is considered as out.
318 if(std::find_if(ptToFill+newNb*nbOfComp,ptToFill+(newNb+1)*nbOfComp,std::bind2nd(std::not_equal_to<double>(),std::numeric_limits<double>::max()))
319 ==ptToFill+(newNb+1)*nbOfComp)
320 std::copy(ptSrc+i*nbOfComp,ptSrc+(i+1)*nbOfComp,ptToFill+newNb*nbOfComp);
323 std::transform(ptSrc+i*nbOfComp,ptSrc+(i+1)*nbOfComp,ptToFill+newNb*nbOfComp,(double *)tmp,std::minus<double>());
324 std::transform((double *)tmp,((double *)tmp)+nbOfComp,(double *)tmp,std::ptr_fun<double,double>(fabs));
325 //if(!std::equal(ptSrc+i*nbOfComp,ptSrc+(i+1)*nbOfComp,ptToFill+newNb*nbOfComp))
326 if(*std::max_element((double *)tmp,((double *)tmp)+nbOfComp)>eps)
329 std::ostringstream oss;
330 oss << msg << " " << i << " and " << std::find(old2NewPtr,old2NewPtr+i,newNb)-old2NewPtr
331 << " have been merged and " << msg << " field on them are different !";
332 throw INTERP_KERNEL::Exception(oss.str().c_str());
340 void MEDCouplingFieldDiscretization::RenumberEntitiesFromN2OArr(const int *new2OldPtr, int new2OldSz, DataArrayDouble *arr, const char *msg)
342 int nbOfComp=arr->getNumberOfComponents();
343 DataArrayDouble *arrCpy=arr->deepCpy();
344 const double *ptSrc=arrCpy->getConstPointer();
345 arr->reAlloc(new2OldSz);
346 double *ptToFill=arr->getPointer();
347 for(int i=0;i<new2OldSz;i++)
349 int oldNb=new2OldPtr[i];
350 std::copy(ptSrc+oldNb*nbOfComp,ptSrc+(oldNb+1)*nbOfComp,ptToFill+i*nbOfComp);
355 MEDCouplingFieldDiscretization::~MEDCouplingFieldDiscretization()
359 TypeOfField MEDCouplingFieldDiscretizationP0::getEnum() const
364 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationP0::clone() const
366 return new MEDCouplingFieldDiscretizationP0;
369 std::string MEDCouplingFieldDiscretizationP0::getStringRepr() const
371 return std::string(REPR);
374 const char *MEDCouplingFieldDiscretizationP0::getRepr() const
379 bool MEDCouplingFieldDiscretizationP0::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
381 const MEDCouplingFieldDiscretizationP0 *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationP0 *>(other);
384 reason="Spatial discrtization of this is ON_CELLS, which is not the case of other.";
388 int MEDCouplingFieldDiscretizationP0::getNumberOfTuples(const MEDCouplingMesh *mesh) const
390 return mesh->getNumberOfCells();
393 int MEDCouplingFieldDiscretizationP0::getNumberOfMeshPlaces(const MEDCouplingMesh *mesh) const
395 return mesh->getNumberOfCells();
398 DataArrayInt *MEDCouplingFieldDiscretizationP0::getOffsetArr(const MEDCouplingMesh *mesh) const
400 int nbOfTuples=mesh->getNumberOfCells();
401 DataArrayInt *ret=DataArrayInt::New();
402 ret->alloc(nbOfTuples+1,1);
407 void MEDCouplingFieldDiscretizationP0::renumberArraysForCell(const MEDCouplingMesh *mesh, const std::vector<DataArrayDouble *>& arrays,
408 const int *old2NewBg, bool check) throw(INTERP_KERNEL::Exception)
410 const int *array=old2NewBg;
412 array=DataArrayInt::CheckAndPreparePermutation(old2NewBg,old2NewBg+mesh->getNumberOfCells());
413 for(std::vector<DataArrayDouble *>::const_iterator it=arrays.begin();it!=arrays.end();it++)
416 (*it)->renumberInPlace(array);
422 DataArrayDouble *MEDCouplingFieldDiscretizationP0::getLocalizationOfDiscValues(const MEDCouplingMesh *mesh) const
424 return mesh->getBarycenterAndOwner();
427 void MEDCouplingFieldDiscretizationP0::computeMeshRestrictionFromTupleIds(const MEDCouplingMesh *mesh, const int *partBg, const int *partEnd,
428 DataArrayInt *&cellRest)
430 cellRest=DataArrayInt::New();
431 cellRest->alloc((int)std::distance(partBg,partEnd),1);
432 std::copy(partBg,partEnd,cellRest->getPointer());
435 void MEDCouplingFieldDiscretizationP0::checkCompatibilityWithNature(NatureOfField nat) const throw(INTERP_KERNEL::Exception)
439 void MEDCouplingFieldDiscretizationP0::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArrayDouble *da) const throw(INTERP_KERNEL::Exception)
441 if(mesh->getNumberOfCells()!=da->getNumberOfTuples())
443 std::ostringstream message;
444 message << "Field on cells invalid because there are " << mesh->getNumberOfCells();
445 message << " cells in mesh and " << da->getNumberOfTuples() << " tuples in field !";
446 throw INTERP_KERNEL::Exception(message.str().c_str());
450 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationP0::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
452 return mesh->getMeasureField(isAbs);
455 void MEDCouplingFieldDiscretizationP0::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
457 int id=mesh->getCellContainingPoint(loc,_precision);
459 throw INTERP_KERNEL::Exception("Specified point is detected outside of mesh : unable to apply P0::getValueOn !");
460 arr->getTuple(id,res);
463 void MEDCouplingFieldDiscretizationP0::getValueOnPos(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, int i, int j, int k, double *res) const
465 const MEDCouplingCMesh *meshC=dynamic_cast<const MEDCouplingCMesh *>(mesh);
467 throw INTERP_KERNEL::Exception("P0::getValueOnPos is only accessible for structured meshes !");
468 int id=meshC->getCellIdFromPos(i,j,k);
469 arr->getTuple(id,res);
472 DataArrayDouble *MEDCouplingFieldDiscretizationP0::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const
474 std::vector<int> elts,eltsIndex;
475 mesh->getCellsContainingPoints(loc,nbOfPoints,_precision,elts,eltsIndex);
476 int spaceDim=mesh->getSpaceDimension();
477 int nbOfComponents=arr->getNumberOfComponents();
478 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
479 ret->alloc(nbOfPoints,nbOfComponents);
480 double *ptToFill=ret->getPointer();
481 for(int i=0;i<nbOfPoints;i++,ptToFill+=nbOfComponents)
482 if(eltsIndex[i+1]-eltsIndex[i]>=1)
483 arr->getTuple(elts[eltsIndex[i]],ptToFill);
486 std::ostringstream oss; oss << "Point #" << i << " with coordinates : (";
487 std::copy(loc+i*spaceDim,loc+(i+1)*spaceDim,std::ostream_iterator<double>(oss,", "));
488 oss << ") detected outside mesh : unable to apply P0::getValueOnMulti ! ";
489 throw INTERP_KERNEL::Exception(oss.str().c_str());
496 * Nothing to do. It's not a bug.
498 void MEDCouplingFieldDiscretizationP0::renumberValuesOnNodes(double , const int *, int newNbOfNodes, DataArrayDouble *) const
502 void MEDCouplingFieldDiscretizationP0::renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const
504 RenumberEntitiesFromO2NArr(epsOnVals,old2New,newSz,arr,"Cell");
507 void MEDCouplingFieldDiscretizationP0::renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const
509 RenumberEntitiesFromN2OArr(new2old,newSz,arr,"Cell");
513 * This method returns a tuple ids selection from cell ids selection [start;end).
514 * This method is called by MEDCouplingFieldDiscretizationP0::buildSubMeshData to return parameter \b di.
515 * Here for P0 it's very simple !
517 * \return a newly allocated array containing ids to select into the DataArrayDouble of the field.
520 DataArrayInt *MEDCouplingFieldDiscretizationP0::computeTupleIdsToSelectFromCellIds(const MEDCouplingMesh *mesh, const int *startCellIds, const int *endCellIds) const
522 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret=DataArrayInt::New();
523 ret->alloc((int)std::distance(startCellIds,endCellIds),1);
524 std::copy(startCellIds,endCellIds,ret->getPointer());
525 ret->incrRef(); return ret;
529 * This method returns a submesh of 'mesh' instance constituting cell ids contained in array defined as an interval [start;end).
530 * @param di is an array returned that specifies entity ids (here cells ids) in mesh 'mesh' of entity in returned submesh.
531 * Example : The first cell id of returned mesh has the (*di)[0] id in 'mesh'
533 MEDCouplingMesh *MEDCouplingFieldDiscretizationP0::buildSubMeshData(const MEDCouplingMesh *mesh, const int *start, const int *end, DataArrayInt *&di) const
535 MEDCouplingMesh *ret=mesh->buildPart(start,end);
536 di=DataArrayInt::New();
537 di->alloc((int)std::distance(start,end),1);
538 int *pt=di->getPointer();
539 std::copy(start,end,pt);
543 int MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuples(const MEDCouplingMesh *mesh) const
545 return mesh->getNumberOfNodes();
548 int MEDCouplingFieldDiscretizationOnNodes::getNumberOfMeshPlaces(const MEDCouplingMesh *mesh) const
550 return mesh->getNumberOfNodes();
554 * Nothing to do here.
556 void MEDCouplingFieldDiscretizationOnNodes::renumberArraysForCell(const MEDCouplingMesh *, const std::vector<DataArrayDouble *>& arrays,
557 const int *old2NewBg, bool check) throw(INTERP_KERNEL::Exception)
561 DataArrayInt *MEDCouplingFieldDiscretizationOnNodes::getOffsetArr(const MEDCouplingMesh *mesh) const
563 int nbOfTuples=mesh->getNumberOfNodes();
564 DataArrayInt *ret=DataArrayInt::New();
565 ret->alloc(nbOfTuples+1,1);
570 DataArrayDouble *MEDCouplingFieldDiscretizationOnNodes::getLocalizationOfDiscValues(const MEDCouplingMesh *mesh) const
572 return mesh->getCoordinatesAndOwner();
575 void MEDCouplingFieldDiscretizationOnNodes::computeMeshRestrictionFromTupleIds(const MEDCouplingMesh *mesh, const int *partBg, const int *partEnd,
576 DataArrayInt *&cellRest)
578 cellRest=mesh->getCellIdsFullyIncludedInNodeIds(partBg,partEnd);
581 void MEDCouplingFieldDiscretizationOnNodes::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArrayDouble *da) const throw(INTERP_KERNEL::Exception)
583 if(mesh->getNumberOfNodes()!=da->getNumberOfTuples())
585 std::ostringstream message;
586 message << "Field on nodes invalid because there are " << mesh->getNumberOfNodes();
587 message << " nodes in mesh and " << da->getNumberOfTuples() << " tuples in field !";
588 throw INTERP_KERNEL::Exception(message.str().c_str());
593 * This method returns a submesh of 'mesh' instance constituting cell ids contained in array defined as an interval [start;end).
594 * @param di is an array returned that specifies entity ids (here nodes ids) in mesh 'mesh' of entity in returned submesh.
595 * Example : The first node id of returned mesh has the (*di)[0] id in 'mesh'
597 MEDCouplingMesh *MEDCouplingFieldDiscretizationOnNodes::buildSubMeshData(const MEDCouplingMesh *mesh, const int *start, const int *end, DataArrayInt *&di) const
599 MEDCouplingMesh *ret=mesh->buildPartAndReduceNodes(start,end,di);
600 DataArrayInt *di2=di->invertArrayO2N2N2O(ret->getNumberOfNodes());
607 * This method returns a tuple ids selection from cell ids selection [start;end).
608 * This method is called by MEDCouplingFieldDiscretizationP0::buildSubMeshData to return parameter \b di.
609 * Here for P1 only nodes fetched by submesh of mesh[startCellIds:endCellIds) is returned !
611 * \return a newly allocated array containing ids to select into the DataArrayDouble of the field.
614 DataArrayInt *MEDCouplingFieldDiscretizationOnNodes::computeTupleIdsToSelectFromCellIds(const MEDCouplingMesh *mesh, const int *startCellIds, const int *endCellIds) const
617 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP1::computeTupleIdsToSelectFromCellIds : null mesh !");
618 const MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> umesh=mesh->buildUnstructured();
619 MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> umesh2=static_cast<MEDCouplingUMesh *>(umesh->buildPartOfMySelf(startCellIds,endCellIds,true));
620 return umesh2->computeFetchedNodeIds();
623 void MEDCouplingFieldDiscretizationOnNodes::renumberValuesOnNodes(double epsOnVals, const int *old2NewPtr, int newNbOfNodes, DataArrayDouble *arr) const
625 RenumberEntitiesFromO2NArr(epsOnVals,old2NewPtr,newNbOfNodes,arr,"Node");
629 * Nothing to do it's not a bug.
631 void MEDCouplingFieldDiscretizationOnNodes::renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const
636 * Nothing to do it's not a bug.
638 void MEDCouplingFieldDiscretizationOnNodes::renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const
642 void MEDCouplingFieldDiscretizationOnNodes::getValueOnPos(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, int i, int j, int k, double *res) const
644 const MEDCouplingCMesh *meshC=dynamic_cast<const MEDCouplingCMesh *>(mesh);
646 throw INTERP_KERNEL::Exception("OnNodes::getValueOnPos(i,j,k) is only accessible for structured meshes !");
647 int id=meshC->getNodeIdFromPos(i,j,k);
648 arr->getTuple(id,res);
651 TypeOfField MEDCouplingFieldDiscretizationP1::getEnum() const
656 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationP1::clone() const
658 return new MEDCouplingFieldDiscretizationP1;
661 std::string MEDCouplingFieldDiscretizationP1::getStringRepr() const
663 return std::string(REPR);
666 const char *MEDCouplingFieldDiscretizationP1::getRepr() const
671 bool MEDCouplingFieldDiscretizationP1::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
673 const MEDCouplingFieldDiscretizationP1 *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationP1 *>(other);
676 reason="Spatial discrtization of this is ON_NODES, which is not the case of other.";
680 void MEDCouplingFieldDiscretizationP1::checkCompatibilityWithNature(NatureOfField nat) const throw(INTERP_KERNEL::Exception)
682 if(nat!=ConservativeVolumic)
683 throw INTERP_KERNEL::Exception("Invalid nature for P1 field : expected ConservativeVolumic !");
686 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationP1::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
688 return mesh->getMeasureFieldOnNode(isAbs);
691 void MEDCouplingFieldDiscretizationP1::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
693 int id=mesh->getCellContainingPoint(loc,_precision);
695 throw INTERP_KERNEL::Exception("Specified point is detected outside of mesh : unable to apply P1::getValueOn !");
696 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(id);
697 if(type!=INTERP_KERNEL::NORM_SEG2 && type!=INTERP_KERNEL::NORM_TRI3 && type!=INTERP_KERNEL::NORM_TETRA4)
698 throw INTERP_KERNEL::Exception("P1 getValueOn is not specified for not simplex cells !");
699 getValueInCell(mesh,id,arr,loc,res);
703 * This method localizes a point defined by 'loc' in a cell with id 'cellId' into mesh 'mesh'.
704 * The result is put into res expected to be of size at least arr->getNumberOfComponents()
706 void MEDCouplingFieldDiscretizationP1::getValueInCell(const MEDCouplingMesh *mesh, int cellId, const DataArrayDouble *arr, const double *loc, double *res) const
708 std::vector<int> conn;
709 std::vector<double> coo;
710 mesh->getNodeIdsOfCell(cellId,conn);
711 for(std::vector<int>::const_iterator iter=conn.begin();iter!=conn.end();iter++)
712 mesh->getCoordinatesOfNode(*iter,coo);
713 int spaceDim=mesh->getSpaceDimension();
714 std::size_t nbOfNodes=conn.size();
715 std::vector<const double *> vec(nbOfNodes);
716 for(std::size_t i=0;i<nbOfNodes;i++)
717 vec[i]=&coo[i*spaceDim];
718 INTERP_KERNEL::AutoPtr<double> tmp=new double[nbOfNodes];
719 INTERP_KERNEL::barycentric_coords(vec,loc,tmp);
720 int sz=arr->getNumberOfComponents();
721 INTERP_KERNEL::AutoPtr<double> tmp2=new double[sz];
722 std::fill(res,res+sz,0.);
723 for(std::size_t i=0;i<nbOfNodes;i++)
725 arr->getTuple(conn[i],(double *)tmp2);
726 std::transform((double *)tmp2,((double *)tmp2)+sz,(double *)tmp2,std::bind2nd(std::multiplies<double>(),tmp[i]));
727 std::transform(res,res+sz,(double *)tmp2,res,std::plus<double>());
731 DataArrayDouble *MEDCouplingFieldDiscretizationP1::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const
733 std::vector<int> elts,eltsIndex;
734 mesh->getCellsContainingPoints(loc,nbOfPoints,_precision,elts,eltsIndex);
735 int spaceDim=mesh->getSpaceDimension();
736 int nbOfComponents=arr->getNumberOfComponents();
737 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
738 ret->alloc(nbOfPoints,nbOfComponents);
739 double *ptToFill=ret->getPointer();
740 for(int i=0;i<nbOfPoints;i++)
741 if(eltsIndex[i+1]-eltsIndex[i]>=1)
742 getValueInCell(mesh,elts[eltsIndex[i]],arr,loc+i*spaceDim,ptToFill+i*nbOfComponents);
745 std::ostringstream oss; oss << "Point #" << i << " with coordinates : (";
746 std::copy(loc+i*spaceDim,loc+(i+1)*spaceDim,std::ostream_iterator<double>(oss,", "));
747 oss << ") detected outside mesh : unable to apply P1::getValueOnMulti ! ";
748 throw INTERP_KERNEL::Exception(oss.str().c_str());
754 MEDCouplingFieldDiscretizationPerCell::MEDCouplingFieldDiscretizationPerCell():_discr_per_cell(0)
758 MEDCouplingFieldDiscretizationPerCell::~MEDCouplingFieldDiscretizationPerCell()
761 _discr_per_cell->decrRef();
764 MEDCouplingFieldDiscretizationPerCell::MEDCouplingFieldDiscretizationPerCell(const MEDCouplingFieldDiscretizationPerCell& other, const int *startCellIds, const int *endCellIds):_discr_per_cell(0)
766 DataArrayInt *arr=other._discr_per_cell;
769 if(startCellIds==0 && endCellIds==0)
770 _discr_per_cell=arr->deepCpy();
772 _discr_per_cell=arr->selectByTupleIdSafe(startCellIds,endCellIds);
776 void MEDCouplingFieldDiscretizationPerCell::updateTime() const
779 updateTimeWith(*_discr_per_cell);
782 void MEDCouplingFieldDiscretizationPerCell::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArrayDouble *da) const throw(INTERP_KERNEL::Exception)
785 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell has no discretization per cell !");
786 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
787 if(nbOfTuples!=mesh->getNumberOfCells())
788 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell has a discretization per cell but it's not matching the underlying mesh !");
791 bool MEDCouplingFieldDiscretizationPerCell::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
793 const MEDCouplingFieldDiscretizationPerCell *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationPerCell *>(other);
796 reason="Spatial discrtization of this is ON_GAUSS, which is not the case of other.";
799 if(_discr_per_cell==0)
800 return otherC->_discr_per_cell==0;
801 if(otherC->_discr_per_cell==0)
803 bool ret=_discr_per_cell->isEqualIfNotWhy(*otherC->_discr_per_cell,reason);
805 reason.insert(0,"Field discretization per cell DataArrayInt given the discid per cell :");
809 bool MEDCouplingFieldDiscretizationPerCell::isEqualWithoutConsideringStr(const MEDCouplingFieldDiscretization *other, double eps) const
811 const MEDCouplingFieldDiscretizationPerCell *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationPerCell *>(other);
814 if(_discr_per_cell==0)
815 return otherC->_discr_per_cell==0;
816 if(otherC->_discr_per_cell==0)
818 return _discr_per_cell->isEqualWithoutConsideringStr(*otherC->_discr_per_cell);
822 * This method is typically the first step of renumbering. The impact on _discr_per_cell is necessary here.
823 * virtualy by this method.
825 void MEDCouplingFieldDiscretizationPerCell::renumberCells(const int *old2NewBg, bool check) throw(INTERP_KERNEL::Exception)
827 int nbCells=_discr_per_cell->getNumberOfTuples();
828 const int *array=old2NewBg;
830 array=DataArrayInt::CheckAndPreparePermutation(old2NewBg,old2NewBg+nbCells);
832 DataArrayInt *dpc=_discr_per_cell->renumber(array);
833 _discr_per_cell->decrRef();
837 delete [] const_cast<int *>(array);
840 void MEDCouplingFieldDiscretizationPerCell::buildDiscrPerCellIfNecessary(const MEDCouplingMesh *m)
844 _discr_per_cell=DataArrayInt::New();
845 int nbTuples=m->getNumberOfCells();
846 _discr_per_cell->alloc(nbTuples,1);
847 int *ptr=_discr_per_cell->getPointer();
848 std::fill(ptr,ptr+nbTuples,DFT_INVALID_LOCID_VALUE);
852 void MEDCouplingFieldDiscretizationPerCell::checkNoOrphanCells() const throw(INTERP_KERNEL::Exception)
855 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell::checkNoOrphanCells : no discretization defined !");
856 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> test=_discr_per_cell->getIdsEqual(DFT_INVALID_LOCID_VALUE);
857 if(test->getNumberOfTuples()!=0)
858 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell::checkNoOrphanCells : presence of orphan cells !");
861 const DataArrayInt *MEDCouplingFieldDiscretizationPerCell::getArrayOfDiscIds() const
863 return _discr_per_cell;
866 MEDCouplingFieldDiscretizationGauss::MEDCouplingFieldDiscretizationGauss()
870 MEDCouplingFieldDiscretizationGauss::MEDCouplingFieldDiscretizationGauss(const MEDCouplingFieldDiscretizationGauss& other, const int *startCellIds, const int *endCellIds):MEDCouplingFieldDiscretizationPerCell(other,startCellIds,endCellIds),_loc(other._loc)
874 TypeOfField MEDCouplingFieldDiscretizationGauss::getEnum() const
879 bool MEDCouplingFieldDiscretizationGauss::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
881 const MEDCouplingFieldDiscretizationGauss *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationGauss *>(other);
884 reason="Spatial discrtization of this is ON_GAUSS, which is not the case of other.";
887 if(!MEDCouplingFieldDiscretizationPerCell::isEqualIfNotWhy(other,eps,reason))
889 if(_loc.size()!=otherC->_loc.size())
891 reason="Gauss spatial discretization : localization sizes differ";
894 std::size_t sz=_loc.size();
895 for(std::size_t i=0;i<sz;i++)
896 if(!_loc[i].isEqual(otherC->_loc[i],eps))
898 std::ostringstream oss; oss << "Gauss spatial discretization : Localization #" << i << " differ from this to other.";
905 bool MEDCouplingFieldDiscretizationGauss::isEqualWithoutConsideringStr(const MEDCouplingFieldDiscretization *other, double eps) const
907 const MEDCouplingFieldDiscretizationGauss *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationGauss *>(other);
910 if(!MEDCouplingFieldDiscretizationPerCell::isEqualWithoutConsideringStr(other,eps))
912 if(_loc.size()!=otherC->_loc.size())
914 std::size_t sz=_loc.size();
915 for(std::size_t i=0;i<sz;i++)
916 if(!_loc[i].isEqual(otherC->_loc[i],eps))
921 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationGauss::clone() const
923 return new MEDCouplingFieldDiscretizationGauss(*this);
926 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationGauss::clonePart(const int *startCellIds, const int *endCellIds) const
928 return new MEDCouplingFieldDiscretizationGauss(*this,startCellIds,endCellIds);
931 std::string MEDCouplingFieldDiscretizationGauss::getStringRepr() const
933 std::ostringstream oss; oss << REPR << "." << std::endl;
936 if(_discr_per_cell->isAllocated())
938 oss << "Discretization per cell : ";
939 std::copy(_discr_per_cell->begin(),_discr_per_cell->end(),std::ostream_iterator<int>(oss,", "));
943 oss << "Presence of " << _loc.size() << " localizations." << std::endl;
945 for(std::vector<MEDCouplingGaussLocalization>::const_iterator it=_loc.begin();it!=_loc.end();it++,i++)
947 oss << "+++++ Localization #" << i << " +++++" << std::endl;
948 oss << (*it).getStringRepr();
949 oss << "++++++++++" << std::endl;
954 const char *MEDCouplingFieldDiscretizationGauss::getRepr() const
959 int MEDCouplingFieldDiscretizationGauss::getNumberOfTuples(const MEDCouplingMesh *) const
962 const int *dcPtr=_discr_per_cell->getConstPointer();
963 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
964 for(const int *w=dcPtr;w!=dcPtr+nbOfTuples;w++)
965 ret+=_loc[*w].getNumberOfGaussPt();
969 int MEDCouplingFieldDiscretizationGauss::getNumberOfMeshPlaces(const MEDCouplingMesh *mesh) const
971 return mesh->getNumberOfCells();
974 DataArrayInt *MEDCouplingFieldDiscretizationGauss::getOffsetArr(const MEDCouplingMesh *mesh) const
976 int nbOfTuples=mesh->getNumberOfCells();
977 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret=DataArrayInt::New();
978 ret->alloc(nbOfTuples+1,1);
979 int *retPtr=ret->getPointer();
980 const int *start=_discr_per_cell->getConstPointer();
981 if(_discr_per_cell->getNumberOfTuples()!=nbOfTuples)
982 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getOffsetArr : mismatch between the mesh and the discretization ids array length !");
983 int maxPossible=(int)_loc.size();
985 for(int i=0;i<nbOfTuples;i++,start++)
987 if(*start>=0 && *start<maxPossible)
988 retPtr[i+1]=retPtr[i]+_loc[*start].getNumberOfGaussPt();
991 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getOffsetArr : At position #" << i << " the locid = " << *start << " whereas it should be in [0," << maxPossible << ") !";
992 throw INTERP_KERNEL::Exception(oss.str().c_str());
999 void MEDCouplingFieldDiscretizationGauss::renumberArraysForCell(const MEDCouplingMesh *mesh, const std::vector<DataArrayDouble *>& arrays,
1000 const int *old2NewBg, bool check) throw(INTERP_KERNEL::Exception)
1002 const int *array=old2NewBg;
1004 array=DataArrayInt::CheckAndPreparePermutation(old2NewBg,old2NewBg+mesh->getNumberOfCells());
1005 int nbOfCells=_discr_per_cell->getNumberOfTuples();
1006 int nbOfTuples=getNumberOfTuples(0);
1007 const int *dcPtr=_discr_per_cell->getConstPointer();
1008 int *array2=new int[nbOfTuples];//stores the final conversion array old2New to give to arrays in renumberInPlace.
1009 int *array3=new int[nbOfCells];//store for each cell in present dcp array (already renumbered) the offset needed by each cell in new numbering.
1011 for(int i=1;i<nbOfCells;i++)
1012 array3[i]=array3[i-1]+_loc[dcPtr[i-1]].getNumberOfGaussPt();
1014 for(int i=0;i<nbOfCells;i++)
1016 int nbOfGaussPt=_loc[dcPtr[array[i]]].getNumberOfGaussPt();
1017 for(int k=0;k<nbOfGaussPt;k++,j++)
1018 array2[j]=array3[array[i]]+k;
1021 for(std::vector<DataArrayDouble *>::const_iterator it=arrays.begin();it!=arrays.end();it++)
1023 (*it)->renumberInPlace(array2);
1026 delete [] const_cast<int*>(array);
1029 DataArrayDouble *MEDCouplingFieldDiscretizationGauss::getLocalizationOfDiscValues(const MEDCouplingMesh *mesh) const
1031 checkNoOrphanCells();
1032 MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> umesh=mesh->buildUnstructured();//in general do nothing
1033 int nbOfTuples=getNumberOfTuples(mesh);
1034 DataArrayDouble *ret=DataArrayDouble::New();
1035 int spaceDim=mesh->getSpaceDimension();
1036 ret->alloc(nbOfTuples,spaceDim);
1037 std::vector< int > locIds;
1038 std::vector<DataArrayInt *> parts=splitIntoSingleGaussDicrPerCellType(locIds);
1039 std::vector< MEDCouplingAutoRefCountObjectPtr<DataArrayInt> > parts2(parts.size());
1040 std::copy(parts.begin(),parts.end(),parts2.begin());
1041 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> offsets=buildNbOfGaussPointPerCellField();
1042 offsets->computeOffsets();
1043 const int *ptrOffsets=offsets->getConstPointer();
1044 const double *coords=umesh->getCoords()->getConstPointer();
1045 const int *connI=umesh->getNodalConnectivityIndex()->getConstPointer();
1046 const int *conn=umesh->getNodalConnectivity()->getConstPointer();
1047 double *valsToFill=ret->getPointer();
1048 for(std::size_t i=0;i<parts2.size();i++)
1050 INTERP_KERNEL::GaussCoords calculator;
1052 const MEDCouplingGaussLocalization& cli=_loc[locIds[i]];//curLocInfo
1053 INTERP_KERNEL::NormalizedCellType typ=cli.getType();
1054 const std::vector<double>& wg=cli.getWeights();
1055 calculator.addGaussInfo(typ,INTERP_KERNEL::CellModel::GetCellModel(typ).getDimension(),
1056 &cli.getGaussCoords()[0],(int)wg.size(),&cli.getRefCoords()[0],
1057 INTERP_KERNEL::CellModel::GetCellModel(typ).getNumberOfNodes());
1059 int nbt=parts2[i]->getNumberOfTuples();
1060 for(const int *w=parts2[i]->getConstPointer();w!=parts2[i]->getConstPointer()+nbt;w++)
1061 calculator.calculateCoords(cli.getType(),coords,spaceDim,conn+connI[*w]+1,valsToFill+spaceDim*(ptrOffsets[*w]));
1063 ret->copyStringInfoFrom(*umesh->getCoords());
1067 void MEDCouplingFieldDiscretizationGauss::computeMeshRestrictionFromTupleIds(const MEDCouplingMesh *mesh, const int *partBg, const int *partEnd,
1068 DataArrayInt *&cellRest)
1070 throw INTERP_KERNEL::Exception("Not implemented yet !");
1076 void MEDCouplingFieldDiscretizationGauss::checkCompatibilityWithNature(NatureOfField nat) const throw(INTERP_KERNEL::Exception)
1080 void MEDCouplingFieldDiscretizationGauss::getTinySerializationIntInformation(std::vector<int>& tinyInfo) const
1084 val=_discr_per_cell->getNumberOfTuples();
1085 tinyInfo.push_back(val);
1086 tinyInfo.push_back((int)_loc.size());
1088 tinyInfo.push_back(-1);
1090 tinyInfo.push_back(_loc[0].getDimension());
1091 for(std::vector<MEDCouplingGaussLocalization>::const_iterator iter=_loc.begin();iter!=_loc.end();iter++)
1092 (*iter).pushTinySerializationIntInfo(tinyInfo);
1095 void MEDCouplingFieldDiscretizationGauss::getTinySerializationDbleInformation(std::vector<double>& tinyInfo) const
1097 for(std::vector<MEDCouplingGaussLocalization>::const_iterator iter=_loc.begin();iter!=_loc.end();iter++)
1098 (*iter).pushTinySerializationDblInfo(tinyInfo);
1101 void MEDCouplingFieldDiscretizationGauss::getSerializationIntArray(DataArrayInt *& arr) const
1105 arr=_discr_per_cell;
1108 void MEDCouplingFieldDiscretizationGauss::resizeForUnserialization(const std::vector<int>& tinyInfo, DataArrayInt *& arr)
1110 int val=tinyInfo[0];
1113 _discr_per_cell=DataArrayInt::New();
1114 _discr_per_cell->alloc(val,1);
1118 arr=_discr_per_cell;
1119 int nbOfLoc=tinyInfo[1];
1121 int dim=tinyInfo[2];
1124 delta=((int)tinyInfo.size()-3)/nbOfLoc;
1125 for(int i=0;i<nbOfLoc;i++)
1127 std::vector<int> tmp(tinyInfo.begin()+3+i*delta,tinyInfo.begin()+3+(i+1)*delta);
1128 MEDCouplingGaussLocalization elt=MEDCouplingGaussLocalization::BuildNewInstanceFromTinyInfo(dim,tmp);
1129 _loc.push_back(elt);
1133 void MEDCouplingFieldDiscretizationGauss::finishUnserialization(const std::vector<double>& tinyInfo)
1135 double *tmp=new double[tinyInfo.size()];
1136 std::copy(tinyInfo.begin(),tinyInfo.end(),tmp);
1137 const double *work=tmp;
1138 for(std::vector<MEDCouplingGaussLocalization>::iterator iter=_loc.begin();iter!=_loc.end();iter++)
1139 work=(*iter).fillWithValues(work);
1143 double MEDCouplingFieldDiscretizationGauss::getIJK(const MEDCouplingMesh *mesh, const DataArrayDouble *da,
1144 int cellId, int nodeIdInCell, int compoId) const throw(INTERP_KERNEL::Exception)
1146 int offset=getOffsetOfCell(cellId);
1147 return da->getIJ(offset+nodeIdInCell,compoId);
1150 void MEDCouplingFieldDiscretizationGauss::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArrayDouble *da) const throw(INTERP_KERNEL::Exception)
1152 MEDCouplingFieldDiscretizationPerCell::checkCoherencyBetween(mesh,da);
1153 for(std::vector<MEDCouplingGaussLocalization>::const_iterator iter=_loc.begin();iter!=_loc.end();iter++)
1154 (*iter).checkCoherency();
1155 int nbOfDesc=(int)_loc.size();
1156 int nbOfCells=mesh->getNumberOfCells();
1157 const int *dc=_discr_per_cell->getConstPointer();
1158 for(int i=0;i<nbOfCells;i++)
1162 std::ostringstream oss; oss << "Cell # " << i << " of mesh \"" << mesh->getName() << "\" has an undefined gauss location ! Should never happend !";
1163 throw INTERP_KERNEL::Exception(oss.str().c_str());
1167 std::ostringstream oss; oss << "Cell # " << i << " of mesh \"" << mesh->getName() << "\" has no gauss location !";
1168 throw INTERP_KERNEL::Exception(oss.str().c_str());
1170 if(mesh->getTypeOfCell(i)!=_loc[dc[i]].getType())
1172 std::ostringstream oss; oss << "Types of mesh and gauss location mismatch for cell # " << i;
1173 throw INTERP_KERNEL::Exception(oss.str().c_str());
1176 int nbOfTuples=getNumberOfTuples(mesh);
1177 if(nbOfTuples!=da->getNumberOfTuples())
1179 std::ostringstream oss; oss << "Invalid number of tuples in the array : expecting " << nbOfTuples << " !";
1180 throw INTERP_KERNEL::Exception(oss.str().c_str());
1184 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationGauss::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
1186 throw INTERP_KERNEL::Exception("Not implemented yet !");
1189 void MEDCouplingFieldDiscretizationGauss::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
1191 throw INTERP_KERNEL::Exception("Not implemented yet !");
1194 void MEDCouplingFieldDiscretizationGauss::getValueOnPos(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, int i, int j, int k, double *res) const
1196 throw INTERP_KERNEL::Exception("getValueOnPos(i,j,k) : Not applyable for Gauss points !");
1199 DataArrayDouble *MEDCouplingFieldDiscretizationGauss::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const
1201 throw INTERP_KERNEL::Exception("getValueOnMulti : Not implemented yet for gauss points !");
1204 MEDCouplingMesh *MEDCouplingFieldDiscretizationGauss::buildSubMeshData(const MEDCouplingMesh *mesh, const int *start, const int *end, DataArrayInt *&di) const
1206 di=computeTupleIdsToSelectFromCellIds(mesh,start,end);
1207 return mesh->buildPart(start,end);
1211 * This method returns a tuple ids selection from cell ids selection [start;end).
1212 * This method is called by MEDCouplingFieldDiscretizationGauss::buildSubMeshData to return parameter \b di.
1214 * \return a newly allocated array containing ids to select into the DataArrayDouble of the field.
1217 DataArrayInt *MEDCouplingFieldDiscretizationGauss::computeTupleIdsToSelectFromCellIds(const MEDCouplingMesh *mesh, const int *startCellIds, const int *endCellIds) const
1220 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::computeTupleIdsToSelectFromCellIds : null mesh !");
1221 if(!_discr_per_cell)
1222 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::computeTupleIdsToSelectFromCellIds : null discretization ids !");
1223 int nbOfCells=mesh->getNumberOfCells();
1224 if(_discr_per_cell->getNumberOfTuples()!=nbOfCells)
1225 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::computeTupleIdsToSelectFromCellIds : mismatch of nb of tuples of cell ids array and number of cells !");
1226 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=DataArrayInt::New(); nbOfNodesPerCell->alloc(nbOfCells,1);
1227 int *retPtr=nbOfNodesPerCell->getPointer();
1228 const int *pt=_discr_per_cell->getConstPointer();
1229 int nbMaxOfLocId=(int)_loc.size();
1230 for(int i=0;i<nbOfCells;i++,retPtr++,pt++)
1232 if(*pt>=0 && *pt<nbMaxOfLocId)
1233 *retPtr=_loc[*pt].getNumberOfGaussPt();
1235 nbOfNodesPerCell->computeOffsets2();
1236 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> sel=DataArrayInt::New(); sel->useArray(startCellIds,false,CPP_DEALLOC,(int)std::distance(startCellIds,endCellIds),1);
1237 return sel->buildExplicitArrByRanges(nbOfNodesPerCell);
1241 * No implementation needed !
1243 void MEDCouplingFieldDiscretizationGauss::renumberValuesOnNodes(double , const int *, int newNbOfNodes, DataArrayDouble *) const
1247 void MEDCouplingFieldDiscretizationGauss::renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const
1249 throw INTERP_KERNEL::Exception("Not implemented yet !");
1252 void MEDCouplingFieldDiscretizationGauss::renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const
1254 throw INTERP_KERNEL::Exception("Number of cells has changed and becomes higher with some cells that have been split ! Unable to conserve the Gauss field !");
1257 void MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnType(const MEDCouplingMesh *m, INTERP_KERNEL::NormalizedCellType type, const std::vector<double>& refCoo,
1258 const std::vector<double>& gsCoo, const std::vector<double>& wg) throw(INTERP_KERNEL::Exception)
1260 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
1261 if((int)cm.getDimension()!=m->getMeshDimension())
1263 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnType : mismatch of dimensions ! MeshDim==" << m->getMeshDimension();
1264 oss << " whereas Type '" << cm.getRepr() << "' has dimension " << cm.getDimension() << " !";
1265 throw INTERP_KERNEL::Exception(oss.str().c_str());
1267 buildDiscrPerCellIfNecessary(m);
1268 int id=(int)_loc.size();
1269 MEDCouplingGaussLocalization elt(type,refCoo,gsCoo,wg);
1270 _loc.push_back(elt);
1271 int *ptr=_discr_per_cell->getPointer();
1272 int nbCells=m->getNumberOfCells();
1273 for(int i=0;i<nbCells;i++)
1274 if(m->getTypeOfCell(i)==type)
1276 zipGaussLocalizations();
1279 void MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnCells(const MEDCouplingMesh *m, const int *begin, const int *end, const std::vector<double>& refCoo,
1280 const std::vector<double>& gsCoo, const std::vector<double>& wg) throw(INTERP_KERNEL::Exception)
1282 buildDiscrPerCellIfNecessary(m);
1283 if(std::distance(begin,end)<1)
1284 throw INTERP_KERNEL::Exception("Size of [begin,end) must be equal or greater than 1 !");
1285 INTERP_KERNEL::NormalizedCellType type=m->getTypeOfCell(*begin);
1286 MEDCouplingGaussLocalization elt(type,refCoo,gsCoo,wg);
1287 int id=(int)_loc.size();
1288 int *ptr=_discr_per_cell->getPointer();
1289 for(const int *w=begin+1;w!=end;w++)
1291 if(m->getTypeOfCell(*w)!=type)
1293 std::ostringstream oss; oss << "The cell with id " << *w << " has been detected to be incompatible in the [begin,end) array specified !";
1294 throw INTERP_KERNEL::Exception(oss.str().c_str());
1298 for(const int *w2=begin;w2!=end;w2++)
1301 _loc.push_back(elt);
1302 zipGaussLocalizations();
1305 void MEDCouplingFieldDiscretizationGauss::clearGaussLocalizations() throw(INTERP_KERNEL::Exception)
1309 _discr_per_cell->decrRef();
1315 MEDCouplingGaussLocalization& MEDCouplingFieldDiscretizationGauss::getGaussLocalization(int locId) throw(INTERP_KERNEL::Exception)
1317 checkLocalizationId(locId);
1321 int MEDCouplingFieldDiscretizationGauss::getNbOfGaussLocalization() const throw(INTERP_KERNEL::Exception)
1323 return (int)_loc.size();
1326 int MEDCouplingFieldDiscretizationGauss::getGaussLocalizationIdOfOneCell(int cellId) const throw(INTERP_KERNEL::Exception)
1328 if(!_discr_per_cell)
1329 throw INTERP_KERNEL::Exception("No Gauss localization still set !");
1330 int locId=_discr_per_cell->getConstPointer()[cellId];
1332 throw INTERP_KERNEL::Exception("No Gauss localization set for the specified cell !");
1336 int MEDCouplingFieldDiscretizationGauss::getGaussLocalizationIdOfOneType(INTERP_KERNEL::NormalizedCellType type) const throw(INTERP_KERNEL::Exception)
1338 std::set<int> ret=getGaussLocalizationIdsOfOneType(type);
1340 throw INTERP_KERNEL::Exception("No gauss discretization found for the specified type !");
1342 throw INTERP_KERNEL::Exception("Several gauss discretizations have been found for the specified type !");
1343 return *ret.begin();
1346 std::set<int> MEDCouplingFieldDiscretizationGauss::getGaussLocalizationIdsOfOneType(INTERP_KERNEL::NormalizedCellType type) const throw(INTERP_KERNEL::Exception)
1348 if(!_discr_per_cell)
1349 throw INTERP_KERNEL::Exception("No Gauss localization still set !");
1352 for(std::vector<MEDCouplingGaussLocalization>::const_iterator iter=_loc.begin();iter!=_loc.end();iter++,id++)
1353 if((*iter).getType()==type)
1358 void MEDCouplingFieldDiscretizationGauss::getCellIdsHavingGaussLocalization(int locId, std::vector<int>& cellIds) const throw(INTERP_KERNEL::Exception)
1360 if(locId<0 || locId>=(int)_loc.size())
1361 throw INTERP_KERNEL::Exception("Invalid locId given : must be in range [0:getNbOfGaussLocalization()) !");
1362 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1363 const int *ptr=_discr_per_cell->getConstPointer();
1364 for(int i=0;i<nbOfTuples;i++)
1366 cellIds.push_back(i);
1369 const MEDCouplingGaussLocalization& MEDCouplingFieldDiscretizationGauss::getGaussLocalization(int locId) const throw(INTERP_KERNEL::Exception)
1371 checkLocalizationId(locId);
1375 void MEDCouplingFieldDiscretizationGauss::checkLocalizationId(int locId) const throw(INTERP_KERNEL::Exception)
1377 if(locId<0 || locId>=(int)_loc.size())
1378 throw INTERP_KERNEL::Exception("Invalid locId given : must be in range [0:getNbOfGaussLocalization()) !");
1381 int MEDCouplingFieldDiscretizationGauss::getOffsetOfCell(int cellId) const throw(INTERP_KERNEL::Exception)
1384 const int *start=_discr_per_cell->getConstPointer();
1385 for(const int *w=start;w!=start+cellId;w++)
1386 ret+=_loc[*w].getNumberOfGaussPt();
1391 * This method do the assumption that there is no orphan cell. If there is an exception is thrown.
1392 * This method makes the assumption too that '_discr_per_cell' is defined. If not an exception is thrown.
1393 * This method returns a newly created array with number of tuples equals to '_discr_per_cell->getNumberOfTuples' and number of components equal to 1.
1394 * The i_th tuple in returned array is the number of gauss point if the corresponding cell.
1396 DataArrayInt *MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField() const throw(INTERP_KERNEL::Exception)
1398 if(!_discr_per_cell)
1399 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField : no discretization array set !");
1400 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1401 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret=DataArrayInt::New();
1402 const int *w=_discr_per_cell->getConstPointer();
1403 ret->alloc(nbOfTuples,1);
1404 int *valsToFill=ret->getPointer();
1405 for(int i=0;i<nbOfTuples;i++,w++)
1406 if(*w!=DFT_INVALID_LOCID_VALUE)
1407 valsToFill[i]=_loc[*w].getNumberOfGaussPt();
1409 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField : orphan cell detected !");
1415 * This method makes the assumption that _discr_per_cell is set.
1416 * This method reduces as much as possible number size of _loc.
1417 * This method is usefull when several set on same cells has been done and that some Gauss Localization are no more used.
1419 void MEDCouplingFieldDiscretizationGauss::zipGaussLocalizations()
1421 const int *start=_discr_per_cell->getConstPointer();
1422 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1423 INTERP_KERNEL::AutoPtr<int> tmp=new int[_loc.size()];
1424 std::fill((int *)tmp,(int *)tmp+_loc.size(),-2);
1425 for(const int *w=start;w!=start+nbOfTuples;w++)
1429 for(int i=0;i<(int)_loc.size();i++)
1432 if(fid==(int)_loc.size())
1435 int *start2=_discr_per_cell->getPointer();
1436 for(int *w2=start2;w2!=start2+nbOfTuples;w2++)
1439 std::vector<MEDCouplingGaussLocalization> tmpLoc;
1440 for(int i=0;i<(int)_loc.size();i++)
1442 tmpLoc.push_back(_loc[tmp[i]]);
1447 * This method is usefull when 'this' describes a field discretization with several gauss discretization on a \b same cell type.
1448 * For example same NORM_TRI3 cells having 6 gauss points and others with 12 gauss points.
1449 * This method returns 2 arrays with same size : the return value and 'locIds' output parameter.
1450 * For a given i into [0,locIds.size) ret[i] represents the set of cell ids of i_th set an locIds[i] represents the set of discretisation of the set.
1451 * The return vector contains a set of newly created instance to deal with.
1452 * The returned vector represents a \b partition of cells ids with a gauss discretization set.
1454 * If no descretization is set in 'this' and exception will be thrown.
1456 std::vector<DataArrayInt *> MEDCouplingFieldDiscretizationGauss::splitIntoSingleGaussDicrPerCellType(std::vector<int>& locIds) const throw(INTERP_KERNEL::Exception)
1458 if(!_discr_per_cell)
1459 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::splitIntoSingleGaussDicrPerCellType : no descretization set !");
1460 return _discr_per_cell->partitionByDifferentValues(locIds);
1463 MEDCouplingFieldDiscretizationGaussNE::MEDCouplingFieldDiscretizationGaussNE()
1467 TypeOfField MEDCouplingFieldDiscretizationGaussNE::getEnum() const
1472 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationGaussNE::clone() const
1474 return new MEDCouplingFieldDiscretizationGaussNE(*this);
1477 std::string MEDCouplingFieldDiscretizationGaussNE::getStringRepr() const
1479 return std::string(REPR);
1482 const char *MEDCouplingFieldDiscretizationGaussNE::getRepr() const
1487 bool MEDCouplingFieldDiscretizationGaussNE::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
1489 const MEDCouplingFieldDiscretizationGaussNE *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationGaussNE *>(other);
1492 reason="Spatial discrtization of this is ON_GAUSS_NE, which is not the case of other.";
1496 int MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuples(const MEDCouplingMesh *mesh) const
1499 int nbOfCells=mesh->getNumberOfCells();
1500 for(int i=0;i<nbOfCells;i++)
1502 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
1503 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
1505 throw INTERP_KERNEL::Exception("Not implemented yet Gauss node on elements for polygons and polyedrons !");
1506 ret+=cm.getNumberOfNodes();
1511 int MEDCouplingFieldDiscretizationGaussNE::getNumberOfMeshPlaces(const MEDCouplingMesh *mesh) const
1513 return mesh->getNumberOfCells();
1516 DataArrayInt *MEDCouplingFieldDiscretizationGaussNE::getOffsetArr(const MEDCouplingMesh *mesh) const
1518 int nbOfTuples=mesh->getNumberOfCells();
1519 DataArrayInt *ret=DataArrayInt::New();
1520 ret->alloc(nbOfTuples+1,1);
1521 int *retPtr=ret->getPointer();
1523 for(int i=0;i<nbOfTuples;i++)
1525 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
1526 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
1528 throw INTERP_KERNEL::Exception("Not implemented yet Gauss node on elements for polygons and polyedrons !");
1529 retPtr[i+1]=retPtr[i]+cm.getNumberOfNodes();
1534 void MEDCouplingFieldDiscretizationGaussNE::renumberArraysForCell(const MEDCouplingMesh *mesh, const std::vector<DataArrayDouble *>& arrays,
1535 const int *old2NewBg, bool check) throw(INTERP_KERNEL::Exception)
1537 const int *array=old2NewBg;
1539 array=DataArrayInt::CheckAndPreparePermutation(old2NewBg,old2NewBg+mesh->getNumberOfCells());
1540 int nbOfCells=mesh->getNumberOfCells();
1541 int nbOfTuples=getNumberOfTuples(mesh);
1542 int *array2=new int[nbOfTuples];//stores the final conversion array old2New to give to arrays in renumberInPlace.
1543 int *array3=new int[nbOfCells];//store for each cell in after renumbering the offset needed by each cell in new numbering.
1545 for(int i=1;i<nbOfCells;i++)
1547 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell((int)std::distance(array,std::find(array,array+nbOfCells,i-1)));
1548 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
1549 array3[i]=array3[i-1]+cm.getNumberOfNodes();
1552 for(int i=0;i<nbOfCells;i++)
1554 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
1555 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
1556 for(int k=0;k<(int)cm.getNumberOfNodes();k++,j++)
1557 array2[j]=array3[array[i]]+k;
1560 for(std::vector<DataArrayDouble *>::const_iterator it=arrays.begin();it!=arrays.end();it++)
1562 (*it)->renumberInPlace(array2);
1565 delete [] const_cast<int *>(array);
1568 DataArrayDouble *MEDCouplingFieldDiscretizationGaussNE::getLocalizationOfDiscValues(const MEDCouplingMesh *mesh) const
1570 throw INTERP_KERNEL::Exception("Not implemented yet !");
1573 void MEDCouplingFieldDiscretizationGaussNE::computeMeshRestrictionFromTupleIds(const MEDCouplingMesh *mesh, const int *partBg, const int *partEnd,
1574 DataArrayInt *&cellRest)
1576 throw INTERP_KERNEL::Exception("Not implemented yet !");
1579 void MEDCouplingFieldDiscretizationGaussNE::checkCompatibilityWithNature(NatureOfField nat) const throw(INTERP_KERNEL::Exception)
1583 double MEDCouplingFieldDiscretizationGaussNE::getIJK(const MEDCouplingMesh *mesh, const DataArrayDouble *da,
1584 int cellId, int nodeIdInCell, int compoId) const throw(INTERP_KERNEL::Exception)
1587 for(int i=0;i<cellId;i++)
1589 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
1590 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
1591 offset+=cm.getNumberOfNodes();
1593 return da->getIJ(offset+nodeIdInCell,compoId);
1596 void MEDCouplingFieldDiscretizationGaussNE::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArrayDouble *da) const throw(INTERP_KERNEL::Exception)
1598 int nbOfTuples=getNumberOfTuples(mesh);
1599 if(nbOfTuples!=da->getNumberOfTuples())
1601 std::ostringstream oss; oss << "Invalid number of tuples in the array : expecting " << nbOfTuples << " !";
1602 throw INTERP_KERNEL::Exception(oss.str().c_str());
1606 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationGaussNE::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
1608 throw INTERP_KERNEL::Exception("Not implemented yet !");
1611 void MEDCouplingFieldDiscretizationGaussNE::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
1613 throw INTERP_KERNEL::Exception("Not implemented yet !");
1616 void MEDCouplingFieldDiscretizationGaussNE::getValueOnPos(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, int i, int j, int k, double *res) const
1618 throw INTERP_KERNEL::Exception("getValueOnPos(i,j,k) : Not applyable for Gauss points !");
1621 DataArrayDouble *MEDCouplingFieldDiscretizationGaussNE::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const
1623 throw INTERP_KERNEL::Exception("getValueOnMulti : Not implemented for Gauss NE !");
1626 MEDCouplingMesh *MEDCouplingFieldDiscretizationGaussNE::buildSubMeshData(const MEDCouplingMesh *mesh, const int *start, const int *end, DataArrayInt *&di) const
1628 di=computeTupleIdsToSelectFromCellIds(mesh,start,end);
1629 return mesh->buildPart(start,end);
1633 * This method returns a tuple ids selection from cell ids selection [start;end).
1634 * This method is called by MEDCouplingFieldDiscretizationGaussNE::buildSubMeshData to return parameter \b di.
1636 * \return a newly allocated array containing ids to select into the DataArrayDouble of the field.
1639 DataArrayInt *MEDCouplingFieldDiscretizationGaussNE::computeTupleIdsToSelectFromCellIds(const MEDCouplingMesh *mesh, const int *startCellIds, const int *endCellIds) const
1642 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::computeTupleIdsToSelectFromCellIds : null mesh !");
1643 const MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> umesh=mesh->buildUnstructured();
1644 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=umesh->computeNbOfNodesPerCell();
1645 nbOfNodesPerCell->computeOffsets2();
1646 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> sel=DataArrayInt::New(); sel->useArray(startCellIds,false,CPP_DEALLOC,(int)std::distance(startCellIds,endCellIds),1);
1647 return sel->buildExplicitArrByRanges(nbOfNodesPerCell);
1651 * No implementation needed !
1653 void MEDCouplingFieldDiscretizationGaussNE::renumberValuesOnNodes(double , const int *, int newNbOfNodes, DataArrayDouble *) const
1657 void MEDCouplingFieldDiscretizationGaussNE::renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const
1659 throw INTERP_KERNEL::Exception("Not implemented yet !");
1662 void MEDCouplingFieldDiscretizationGaussNE::renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const
1664 throw INTERP_KERNEL::Exception("Not implemented yet !");
1667 MEDCouplingFieldDiscretizationGaussNE::MEDCouplingFieldDiscretizationGaussNE(const MEDCouplingFieldDiscretizationGaussNE& other):MEDCouplingFieldDiscretization(other)
1671 TypeOfField MEDCouplingFieldDiscretizationKriging::getEnum() const
1676 const char *MEDCouplingFieldDiscretizationKriging::getRepr() const
1681 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationKriging::clone() const
1683 return new MEDCouplingFieldDiscretizationKriging;
1686 std::string MEDCouplingFieldDiscretizationKriging::getStringRepr() const
1688 return std::string(REPR);
1691 void MEDCouplingFieldDiscretizationKriging::checkCompatibilityWithNature(NatureOfField nat) const throw(INTERP_KERNEL::Exception)
1693 if(nat!=ConservativeVolumic)
1694 throw INTERP_KERNEL::Exception("Invalid nature for Kriging field : expected ConservativeVolumic !");
1697 bool MEDCouplingFieldDiscretizationKriging::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
1699 const MEDCouplingFieldDiscretizationKriging *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationKriging *>(other);
1702 reason="Spatial discrtization of this is ON_NODES_KR, which is not the case of other.";
1706 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationKriging::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
1708 throw INTERP_KERNEL::Exception("getMeasureField on FieldDiscretizationKriging : not implemented yet !");
1711 void MEDCouplingFieldDiscretizationKriging::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
1713 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> res2=MEDCouplingFieldDiscretizationKriging::getValueOnMulti(arr,mesh,loc,1);
1714 std::copy(res2->begin(),res2->end(),res);
1717 DataArrayDouble *MEDCouplingFieldDiscretizationKriging::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfTargetPoints) const
1719 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coords=getLocalizationOfDiscValues(mesh);
1720 int nbOfPts=coords->getNumberOfTuples();
1721 int dimension=coords->getNumberOfComponents();
1724 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> KnewiK=computeVectorOfCoefficients(mesh,arr,delta);
1726 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> locArr=DataArrayDouble::New();
1727 locArr->useArray(loc,false,CPP_DEALLOC,nbOfTargetPoints,dimension);
1728 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrix2=coords->buildEuclidianDistanceDenseMatrixWith(locArr);
1729 operateOnDenseMatrix(mesh->getSpaceDimension(),nbOfPts*nbOfTargetPoints,matrix2->getPointer());
1730 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrix3=DataArrayDouble::New();
1731 matrix3->alloc((nbOfPts+delta)*nbOfTargetPoints,1);
1732 double *work=matrix3->getPointer();
1733 const double *workCst=matrix2->getConstPointer();
1734 const double *workCst2=loc;
1735 for(int i=0;i<nbOfTargetPoints;i++,workCst+=nbOfPts,workCst2+=delta-1)
1737 for(int j=0;j<nbOfPts;j++)
1738 work[j*nbOfTargetPoints+i]=workCst[j];
1739 work[nbOfPts*nbOfTargetPoints+i]=1.0;
1740 for(int j=0;j<delta-1;j++)
1741 work[(nbOfPts+1+j)*nbOfTargetPoints+i]=workCst2[j];
1744 int nbOfCompo=arr->getNumberOfComponents();
1745 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
1746 ret->alloc(nbOfTargetPoints,nbOfCompo);
1747 INTERP_KERNEL::matrixProduct(KnewiK->getConstPointer(),1,nbOfPts+delta,matrix3->getConstPointer(),nbOfPts+delta,nbOfTargetPoints*nbOfCompo,ret->getPointer());
1753 * This method computes coefficients to apply to each representing points of \a mesh, that is to say the nodes of \a mesh given a field array \a arr whose
1754 * number of tuples should be equal to the number of representing points in \a mesh.
1756 * \param [in] mesh is the sources of nodes on which kriging will be done regarding the parameters and the value of \c this->getSpaceDimension()
1757 * \param [in] arr input field DataArrayDouble whose number of tuples must be equal to the number of nodes in \a mesh
1758 * \param [out] isDrift return if drift coefficients are present in the returned vector of coefficients, and if. If different from 0 there is presence of drift coefficients.
1759 * Whatever the value of \a isDrift the number of tuples of returned DataArrayDouble will be equal to \c arr->getNumberOfTuples() + \a isDrift.
1760 * \return a newly allocated array containing coefficients including or not drift coefficient at the end depending the value of \a isDrift parameter.
1762 DataArrayDouble *MEDCouplingFieldDiscretizationKriging::computeVectorOfCoefficients(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, int& isDrift) const
1764 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coords=getLocalizationOfDiscValues(mesh);
1765 int nbOfPts=coords->getNumberOfTuples();
1766 //int dimension=coords->getNumberOfComponents();
1767 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrix=coords->buildEuclidianDistanceDenseMatrix();
1768 operateOnDenseMatrix(mesh->getSpaceDimension(),nbOfPts*nbOfPts,matrix->getPointer());
1770 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrixWithDrift=performDrift(matrix,coords,isDrift);
1771 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrixInv=DataArrayDouble::New();
1772 matrixInv->alloc((nbOfPts+isDrift)*(nbOfPts+isDrift),1);
1773 INTERP_KERNEL::inverseMatrix(matrixWithDrift->getConstPointer(),nbOfPts+isDrift,matrixInv->getPointer());
1775 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> KnewiK=DataArrayDouble::New();
1776 KnewiK->alloc((nbOfPts+isDrift)*1,1);
1777 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arr2=DataArrayDouble::New();
1778 arr2->alloc((nbOfPts+isDrift)*1,1);
1779 double *work=std::copy(arr->begin(),arr->end(),arr2->getPointer());
1780 std::fill(work,work+isDrift,0.);
1781 INTERP_KERNEL::matrixProduct(matrixInv->getConstPointer(),nbOfPts+isDrift,nbOfPts+isDrift,arr2->getConstPointer(),nbOfPts+isDrift,1,KnewiK->getPointer());
1787 * Apply \f f(x) on each element x in \a matrixPtr. \a matrixPtr is expected to be a dense matrix represented by a chunck of memory of size at least equal to \a nbOfElems.
1789 * \param [in] spaceDimension space dimension of the input mesh on which the Kriging has to be performed
1790 * \param [in] nbOfElems is the result of the product of nb of rows and the nb of columns of matrix \a matrixPtr
1791 * \param [in,out] matrixPtr is the dense matrix whose on each values the operation will be applied
1793 void MEDCouplingFieldDiscretizationKriging::operateOnDenseMatrix(int spaceDimension, int nbOfElems, double *matrixPtr) const
1795 switch(spaceDimension)
1799 for(int i=0;i<nbOfElems;i++)
1801 double val=matrixPtr[i];
1802 matrixPtr[i]=val*val*val;
1807 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationKriging::operateOnDenseMatrix : only dimension 1 implemented !");
1812 * Starting from a square matrix \a matr, this method returns a newly allocated dense square matrix whose \a matr is included in returned matrix
1813 * in the top left corner, and in the remaining returned matrix the parameters to take into account about the kriging drift.
1814 * For the moment only linear srift is implemented.
1816 * \param [in] arr the position of points were input mesh geometry is considered for Kriging
1817 * \param [in] matr input matrix whose drift part will be added
1818 * \param [out] delta the difference between the size of the output matrix and the input matrix \a matr.
1819 * \return a newly allocated matrix bigger than input matrix \a matr.
1821 DataArrayDouble *MEDCouplingFieldDiscretizationKriging::performDrift(const DataArrayDouble *matr, const DataArrayDouble *arr, int& delta) const
1823 int spaceDimension=arr->getNumberOfComponents();
1824 delta=spaceDimension+1;
1825 int szOfMatrix=arr->getNumberOfTuples();
1826 if(szOfMatrix*szOfMatrix!=matr->getNumberOfTuples())
1827 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationKriging::performDrift : invalid size");
1828 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
1829 ret->alloc((szOfMatrix+delta)*(szOfMatrix+delta),1);
1830 const double *srcWork=matr->getConstPointer();
1831 const double *srcWork2=arr->getConstPointer();
1832 double *destWork=ret->getPointer();
1833 for(int i=0;i<szOfMatrix;i++)
1835 destWork=std::copy(srcWork,srcWork+szOfMatrix,destWork);
1836 srcWork+=szOfMatrix;
1838 destWork=std::copy(srcWork2,srcWork2+spaceDimension,destWork);
1839 srcWork2+=spaceDimension;
1841 std::fill(destWork,destWork+szOfMatrix,1.); destWork+=szOfMatrix;
1842 std::fill(destWork,destWork+spaceDimension+1,0.); destWork+=spaceDimension+1;
1843 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arrNoI=arr->toNoInterlace();
1844 srcWork2=arrNoI->getConstPointer();
1845 for(int i=0;i<spaceDimension;i++)
1847 destWork=std::copy(srcWork2,srcWork2+szOfMatrix,destWork);
1848 srcWork2+=szOfMatrix;
1849 std::fill(destWork,destWork+spaceDimension+1,0.);
1850 destWork+=spaceDimension;