1 // Copyright (C) 2007-2013 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"
41 using namespace ParaMEDMEM;
43 const double MEDCouplingFieldDiscretization::DFLT_PRECISION=1.e-12;
45 const char MEDCouplingFieldDiscretizationP0::REPR[]="P0";
47 const TypeOfField MEDCouplingFieldDiscretizationP0::TYPE=ON_CELLS;
49 const char MEDCouplingFieldDiscretizationP1::REPR[]="P1";
51 const TypeOfField MEDCouplingFieldDiscretizationP1::TYPE=ON_NODES;
53 const int MEDCouplingFieldDiscretizationPerCell::DFT_INVALID_LOCID_VALUE=-1;
55 const char MEDCouplingFieldDiscretizationGauss::REPR[]="GAUSS";
57 const TypeOfField MEDCouplingFieldDiscretizationGauss::TYPE=ON_GAUSS_PT;
59 const char MEDCouplingFieldDiscretizationGaussNE::REPR[]="GSSNE";
61 const TypeOfField MEDCouplingFieldDiscretizationGaussNE::TYPE=ON_GAUSS_NE;
63 const char MEDCouplingFieldDiscretizationKriging::REPR[]="KRIGING";
65 const TypeOfField MEDCouplingFieldDiscretizationKriging::TYPE=ON_NODES_KR;
67 // doc is here http://www.code-aster.org/V2/doc/default/fr/man_r/r3/r3.01.01.pdf
68 const double MEDCouplingFieldDiscretizationGaussNE::FGP_SEG2[2]={1.,1.};
69 const double MEDCouplingFieldDiscretizationGaussNE::FGP_SEG3[3]={0.5555555555555556,0.5555555555555556,0.8888888888888888};
70 const double MEDCouplingFieldDiscretizationGaussNE::FGP_SEG4[4]={0.347854845137454,0.347854845137454,0.652145154862546,0.652145154862546};
71 const double MEDCouplingFieldDiscretizationGaussNE::FGP_TRI3[3]={0.16666666666666666,0.16666666666666666,0.16666666666666666};
72 const double MEDCouplingFieldDiscretizationGaussNE::FGP_TRI6[6]={0.0549758718227661,0.0549758718227661,0.0549758718227661,0.11169079483905,0.11169079483905,0.11169079483905};
73 const double MEDCouplingFieldDiscretizationGaussNE::FGP_TRI7[7]={0.062969590272413,0.062969590272413,0.062969590272413,0.066197076394253,0.066197076394253,0.066197076394253,0.1125};
74 const double MEDCouplingFieldDiscretizationGaussNE::FGP_QUAD4[4]={1.,1.,1.,1.};
75 const double MEDCouplingFieldDiscretizationGaussNE::FGP_QUAD9[9]={0.30864197530864196,0.30864197530864196,0.30864197530864196,0.30864197530864196,0.49382716049382713,0.49382716049382713,0.49382716049382713,0.49382716049382713,0.7901234567901234};
76 const double MEDCouplingFieldDiscretizationGaussNE::FGP_TETRA4[4]={0.041666666666666664,0.041666666666666664,0.041666666666666664,0.041666666666666664};
77 const double MEDCouplingFieldDiscretizationGaussNE::FGP_PENTA6[6]={0.16666666666666666,0.16666666666666666,0.16666666666666666,0.16666666666666666,0.16666666666666666,0.16666666666666666};
78 const double MEDCouplingFieldDiscretizationGaussNE::FGP_HEXA8[8]={1.,1.,1.,1.,1.,1.,1.,1.};
79 const double MEDCouplingFieldDiscretizationGaussNE::FGP_HEXA27[27]={0.1714677640603567,0.1714677640603567,0.1714677640603567,0.1714677640603567,0.1714677640603567,0.1714677640603567,0.1714677640603567,0.1714677640603567,0.27434842249657065,0.27434842249657065,0.27434842249657065,0.27434842249657065,0.27434842249657065,0.27434842249657065,0.27434842249657065,0.27434842249657065,0.27434842249657065,0.27434842249657065,0.27434842249657065,0.27434842249657065,0.43895747599451296,0.43895747599451296,0.43895747599451296,0.43895747599451296,0.43895747599451296,0.43895747599451296,0.7023319615912208};
80 const double MEDCouplingFieldDiscretizationGaussNE::FGP_PYRA5[5]={0.13333333333333333,0.13333333333333333,0.13333333333333333,0.13333333333333333,0.13333333333333333};
81 const double MEDCouplingFieldDiscretizationGaussNE::REF_SEG2[2]={-1.,1.};
82 const double MEDCouplingFieldDiscretizationGaussNE::REF_SEG3[3]={-1.,0.,1.};
83 const double MEDCouplingFieldDiscretizationGaussNE::REF_SEG4[4]={-1.,1.,-0.3333333333333333,0.3333333333333333};
84 const double MEDCouplingFieldDiscretizationGaussNE::REF_TRI3[6]={0.,0.,1.,0.,0.,1.};
85 const double MEDCouplingFieldDiscretizationGaussNE::REF_TRI6[12]={0.,0.,1.,0.,0.,1.,0.5,0.,0.5,0.5,0.,0.5};
86 const double MEDCouplingFieldDiscretizationGaussNE::REF_TRI7[14]={0.,0.,1.,0.,0.,1.,0.5,0.,0.5,0.5,0.,0.5,0.3333333333333333,0.3333333333333333};
87 const double MEDCouplingFieldDiscretizationGaussNE::REF_QUAD4[8]={-1.,-1.,1.,-1.,1.,1.,-1.,1.};
88 const double MEDCouplingFieldDiscretizationGaussNE::REF_QUAD8[16]={-1.,-1.,1.,-1.,1.,1.,-1.,1.,0.,-1.,1.,0.,0.,1.,-1.,0.};
89 const double MEDCouplingFieldDiscretizationGaussNE::REF_QUAD9[18]={-1.,-1.,1.,-1.,1.,1.,-1.,1.,0.,-1.,1.,0.,0.,1.,-1.,0.,0.,0.};
90 const double MEDCouplingFieldDiscretizationGaussNE::REF_TETRA4[12]={0.,1.,0.,0.,0.,1.,0.,0.,0.,1.,0.,0.};
91 const double MEDCouplingFieldDiscretizationGaussNE::REF_TETRA10[30]={0.,1.,0.,0.,0.,1.,0.,0.,0.,1.,0.,0.,0.,0.5,0.5,0.,0.,0.5,0.,0.5,0.,0.5,0.5,0.,0.5,0.,0.5,0.5,0.,0.};
92 const double MEDCouplingFieldDiscretizationGaussNE::REF_PENTA6[18]={-1.,1.,0.,-1.,0.,1.,-1.,0.,0.,1.,1.,0.,1.,0.,1.,1.,0.,0.};
93 const double MEDCouplingFieldDiscretizationGaussNE::REF_PENTA15[45]={-1.,1.,0.,-1.,0.,1.,-1.,0.,0.,1.,1.,0.,1.,0.,1.,1.,0.,0.,-1.,0.5,0.5,-1.,0.,0.5,-1.,0.5,0.,0.,1.,0.,0.,0.,1.,0.,0.,0.,1.,0.5,0.5,1.,0.,0.5,1.,0.5,0.};
94 const double MEDCouplingFieldDiscretizationGaussNE::REF_HEXA8[24]={-1.,-1.,-1.,1.,-1.,-1.,1.,1.,-1.,-1.,1.,-1.,-1.,-1.,1.,1.,-1.,1.,1.,1.,1.,-1.,1.,1.};
95 const double MEDCouplingFieldDiscretizationGaussNE::REF_HEXA20[60]={-1.,-1.,-1.,1.,-1.,-1.,1.,1.,-1.,-1.,1.,-1.,-1.,-1.,1.,1.,-1.,1.,1.,1.,1.,-1.,1.,1.,0.,-1.,-1.,1.,0.,-1.,0.,1.,-1.,-1.,0.,-1.,-1.,-1.,0.,1.,-1.,0.,1.,1.,0.,-1.,1.,0.,0.,-1.,1.,1.,0.,1.,0.,1.,1.,-1.,0.,1.};
96 const double MEDCouplingFieldDiscretizationGaussNE::REF_HEXA27[81]={-1.,-1.,-1.,1.,-1.,-1.,1.,1.,-1.,-1.,1.,-1.,-1.,-1.,1.,1.,-1.,1.,1.,1.,1.,-1.,1.,1.,0.,-1.,-1.,1.,0.,-1.,0.,1.,-1.,-1.,0.,-1.,-1.,-1.,0.,1.,-1.,0.,1.,1.,0.,-1.,1.,0.,0.,-1.,1.,1.,0.,1.,0.,1.,1.,-1.,0.,1.,0.,0.,-1.,0.,-1.,0.,1.,0.,0.,0.,1.,0.,-1.,0.,0.,0.,0.,1.,0.,0.,0.};
97 const double MEDCouplingFieldDiscretizationGaussNE::REF_PYRA5[15]={1.,0.,0.,0.,1.,0.,-1.,0.,0.,0.,-1.,0.,0.,0.,1.};
98 const double MEDCouplingFieldDiscretizationGaussNE::REF_PYRA13[39]={1.,0.,0.,0.,1.,0.,-1.,0.,0.,0.,-1.,0.,0.,0.,1.,0.5,0.5,0.,-0.5,0.5,0.,-0.5,-0.5,0.,0.5,-0.5,0.,0.5,0.,0.5,0.,0.5,0.5,-0.5,0.,0.5,0.,-0.5,0.5};
100 MEDCouplingFieldDiscretization::MEDCouplingFieldDiscretization():_precision(DFLT_PRECISION)
104 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretization::New(TypeOfField type)
108 case MEDCouplingFieldDiscretizationP0::TYPE:
109 return new MEDCouplingFieldDiscretizationP0;
110 case MEDCouplingFieldDiscretizationP1::TYPE:
111 return new MEDCouplingFieldDiscretizationP1;
112 case MEDCouplingFieldDiscretizationGauss::TYPE:
113 return new MEDCouplingFieldDiscretizationGauss;
114 case MEDCouplingFieldDiscretizationGaussNE::TYPE:
115 return new MEDCouplingFieldDiscretizationGaussNE;
116 case MEDCouplingFieldDiscretizationKriging::TYPE:
117 return new MEDCouplingFieldDiscretizationKriging;
119 throw INTERP_KERNEL::Exception("Choosen discretization is not implemented yet.");
123 TypeOfField MEDCouplingFieldDiscretization::GetTypeOfFieldFromStringRepr(const char *repr) throw(INTERP_KERNEL::Exception)
125 std::string reprCpp(repr);
126 if(reprCpp==MEDCouplingFieldDiscretizationP0::REPR)
127 return MEDCouplingFieldDiscretizationP0::TYPE;
128 if(reprCpp==MEDCouplingFieldDiscretizationP1::REPR)
129 return MEDCouplingFieldDiscretizationP1::TYPE;
130 if(reprCpp==MEDCouplingFieldDiscretizationGauss::REPR)
131 return MEDCouplingFieldDiscretizationGauss::TYPE;
132 if(reprCpp==MEDCouplingFieldDiscretizationGaussNE::REPR)
133 return MEDCouplingFieldDiscretizationGaussNE::TYPE;
134 if(reprCpp==MEDCouplingFieldDiscretizationKriging::REPR)
135 return MEDCouplingFieldDiscretizationKriging::TYPE;
136 throw INTERP_KERNEL::Exception("Representation does not match with any field discretization !");
139 bool MEDCouplingFieldDiscretization::isEqual(const MEDCouplingFieldDiscretization *other, double eps) const
142 return isEqualIfNotWhy(other,eps,reason);
145 bool MEDCouplingFieldDiscretization::isEqualWithoutConsideringStr(const MEDCouplingFieldDiscretization *other, double eps) const
147 return isEqual(other,eps);
151 * This method is an alias of MEDCouplingFieldDiscretization::clone. It is only here for coherency with all the remaining of MEDCoupling.
152 * \sa MEDCouplingFieldDiscretization::clone.
154 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretization::deepCpy() const
160 * For all field discretization excepted GaussPts the [ \a startCellIds, \a endCellIds ) has no impact on the cloned instance.
162 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretization::clonePart(const int *startCellIds, const int *endCellIds) const
168 * For all field discretization excepted GaussPts the slice( \a beginCellId, \a endCellIds, \a stepCellId ) has no impact on the cloned instance.
170 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretization::clonePartRange(int beginCellIds, int endCellIds, int stepCellIds) const
176 * Excepted for MEDCouplingFieldDiscretizationPerCell no underlying TimeLabel object : nothing to do in generally.
178 void MEDCouplingFieldDiscretization::updateTime() const
182 std::size_t MEDCouplingFieldDiscretization::getHeapMemorySize() const
188 * Computes normL1 of DataArrayDouble instance arr.
189 * @param res output parameter expected to be of size arr->getNumberOfComponents();
190 * @throw when the field discretization fails on getMeasure fields (gauss points for example)
192 void MEDCouplingFieldDiscretization::normL1(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, double *res) const throw(INTERP_KERNEL::Exception)
194 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=getMeasureField(mesh,true);
195 int nbOfCompo=arr->getNumberOfComponents();
196 int nbOfElems=getNumberOfTuples(mesh);
197 std::fill(res,res+nbOfCompo,0.);
198 const double *arrPtr=arr->getConstPointer();
199 const double *volPtr=vol->getArray()->getConstPointer();
201 for(int i=0;i<nbOfElems;i++)
203 double v=fabs(volPtr[i]);
204 for(int j=0;j<nbOfCompo;j++)
205 res[j]+=fabs(arrPtr[i*nbOfCompo+j])*v;
208 std::transform(res,res+nbOfCompo,res,std::bind2nd(std::multiplies<double>(),1./deno));
212 * Computes normL2 of DataArrayDouble instance arr.
213 * @param res output parameter expected to be of size arr->getNumberOfComponents();
214 * @throw when the field discretization fails on getMeasure fields (gauss points for example)
216 void MEDCouplingFieldDiscretization::normL2(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, double *res) const throw(INTERP_KERNEL::Exception)
218 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=getMeasureField(mesh,true);
219 int nbOfCompo=arr->getNumberOfComponents();
220 int nbOfElems=getNumberOfTuples(mesh);
221 std::fill(res,res+nbOfCompo,0.);
222 const double *arrPtr=arr->getConstPointer();
223 const double *volPtr=vol->getArray()->getConstPointer();
225 for(int i=0;i<nbOfElems;i++)
227 double v=fabs(volPtr[i]);
228 for(int j=0;j<nbOfCompo;j++)
229 res[j]+=arrPtr[i*nbOfCompo+j]*arrPtr[i*nbOfCompo+j]*v;
232 std::transform(res,res+nbOfCompo,res,std::bind2nd(std::multiplies<double>(),1./deno));
233 std::transform(res,res+nbOfCompo,res,std::ptr_fun<double,double>(std::sqrt));
237 * Computes integral of DataArrayDouble instance arr.
238 * @param res output parameter expected to be of size arr->getNumberOfComponents();
239 * @throw when the field discretization fails on getMeasure fields (gauss points for example)
241 void MEDCouplingFieldDiscretization::integral(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, bool isWAbs, double *res) const throw(INTERP_KERNEL::Exception)
244 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretization::integral : mesh is NULL !");
246 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretization::integral : input array is NULL !");
247 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=getMeasureField(mesh,isWAbs);
248 int nbOfCompo=arr->getNumberOfComponents();
249 int nbOfElems=getNumberOfTuples(mesh);
250 if(nbOfElems!=arr->getNumberOfTuples())
252 std::ostringstream oss; oss << "MEDCouplingFieldDiscretization::integral : field is not correct ! number of tuples in array is " << arr->getNumberOfTuples();
253 oss << " whereas number of tuples expected is " << nbOfElems << " !";
254 throw INTERP_KERNEL::Exception(oss.str().c_str());
256 std::fill(res,res+nbOfCompo,0.);
257 const double *arrPtr=arr->getConstPointer();
258 const double *volPtr=vol->getArray()->getConstPointer();
259 INTERP_KERNEL::AutoPtr<double> tmp=new double[nbOfCompo];
260 for (int i=0;i<nbOfElems;i++)
262 std::transform(arrPtr+i*nbOfCompo,arrPtr+(i+1)*nbOfCompo,(double *)tmp,std::bind2nd(std::multiplies<double>(),volPtr[i]));
263 std::transform((double *)tmp,(double *)tmp+nbOfCompo,res,res,std::plus<double>());
268 * This method is strictly equivalent to MEDCouplingFieldDiscretization::buildSubMeshData except that it is optimized for input defined as a range of cell ids.
270 * \param [out] beginOut Valid only if \a di is NULL
271 * \param [out] endOut Valid only if \a di is NULL
272 * \param [out] stepOut Valid only if \a di is NULL
273 * \param [out] di is an array returned that specifies entity ids (nodes, cells, Gauss points... ) in array if no output range is foundable.
275 * \sa MEDCouplingFieldDiscretization::buildSubMeshData
277 MEDCouplingMesh *MEDCouplingFieldDiscretization::buildSubMeshDataRange(const MEDCouplingMesh *mesh, int beginCellIds, int endCellIds, int stepCellIds, int& beginOut, int& endOut, int& stepOut, DataArrayInt *&di) const
279 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> da=DataArrayInt::Range(beginCellIds,endCellIds,stepCellIds);
280 return buildSubMeshData(mesh,da->begin(),da->end(),di);
283 void MEDCouplingFieldDiscretization::getSerializationIntArray(DataArrayInt *& arr) const
291 void MEDCouplingFieldDiscretization::getTinySerializationIntInformation(std::vector<int>& tinyInfo) const
298 void MEDCouplingFieldDiscretization::getTinySerializationDbleInformation(std::vector<double>& tinyInfo) const
302 void MEDCouplingFieldDiscretization::resizeForUnserialization(const std::vector<int>& tinyInfo, DataArrayInt *& arr)
310 void MEDCouplingFieldDiscretization::finishUnserialization(const std::vector<double>& tinyInfo)
315 * This method is typically the first step of renumbering. The implementation is empty it is not a bug only gauss is impacted
316 * virtualy by this method.
318 void MEDCouplingFieldDiscretization::renumberCells(const int *old2NewBg, bool check) throw(INTERP_KERNEL::Exception)
322 double MEDCouplingFieldDiscretization::getIJK(const MEDCouplingMesh *mesh, const DataArrayDouble *da,
323 int cellId, int nodeIdInCell, int compoId) const throw(INTERP_KERNEL::Exception)
325 throw INTERP_KERNEL::Exception("getIJK Invalid ! only for GaussPoint and GaussNE discretizations !");
328 void MEDCouplingFieldDiscretization::setGaussLocalizationOnType(const MEDCouplingMesh *m, INTERP_KERNEL::NormalizedCellType type, const std::vector<double>& refCoo,
329 const std::vector<double>& gsCoo, const std::vector<double>& wg) throw(INTERP_KERNEL::Exception)
331 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
334 void MEDCouplingFieldDiscretization::setGaussLocalizationOnCells(const MEDCouplingMesh *m, const int *begin, const int *end, const std::vector<double>& refCoo,
335 const std::vector<double>& gsCoo, const std::vector<double>& wg) throw(INTERP_KERNEL::Exception)
337 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
340 void MEDCouplingFieldDiscretization::clearGaussLocalizations() throw(INTERP_KERNEL::Exception)
342 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
345 MEDCouplingGaussLocalization& MEDCouplingFieldDiscretization::getGaussLocalization(int locId) throw(INTERP_KERNEL::Exception)
347 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
350 const MEDCouplingGaussLocalization& MEDCouplingFieldDiscretization::getGaussLocalization(int locId) const throw(INTERP_KERNEL::Exception)
352 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
355 int MEDCouplingFieldDiscretization::getNbOfGaussLocalization() const throw(INTERP_KERNEL::Exception)
357 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
360 int MEDCouplingFieldDiscretization::getGaussLocalizationIdOfOneCell(int cellId) const throw(INTERP_KERNEL::Exception)
362 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
365 int MEDCouplingFieldDiscretization::getGaussLocalizationIdOfOneType(INTERP_KERNEL::NormalizedCellType type) const throw(INTERP_KERNEL::Exception)
367 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
370 std::set<int> MEDCouplingFieldDiscretization::getGaussLocalizationIdsOfOneType(INTERP_KERNEL::NormalizedCellType type) const throw(INTERP_KERNEL::Exception)
372 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
375 void MEDCouplingFieldDiscretization::getCellIdsHavingGaussLocalization(int locId, std::vector<int>& cellIds) const throw(INTERP_KERNEL::Exception)
377 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
380 void MEDCouplingFieldDiscretization::RenumberEntitiesFromO2NArr(double eps, const int *old2NewPtr, int newNbOfEntity, DataArrayDouble *arr, const char *msg)
383 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretization::RenumberEntitiesFromO2NArr : input array is NULL !");
384 int oldNbOfElems=arr->getNumberOfTuples();
385 int nbOfComp=arr->getNumberOfComponents();
386 int newNbOfTuples=newNbOfEntity;
387 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arrCpy=arr->deepCpy();
388 const double *ptSrc=arrCpy->getConstPointer();
389 arr->reAlloc(newNbOfTuples);
390 double *ptToFill=arr->getPointer();
391 std::fill(ptToFill,ptToFill+nbOfComp*newNbOfTuples,std::numeric_limits<double>::max());
392 INTERP_KERNEL::AutoPtr<double> tmp=new double[nbOfComp];
393 for(int i=0;i<oldNbOfElems;i++)
395 int newNb=old2NewPtr[i];
396 if(newNb>=0)//if newNb<0 the node is considered as out.
398 if(std::find_if(ptToFill+newNb*nbOfComp,ptToFill+(newNb+1)*nbOfComp,std::bind2nd(std::not_equal_to<double>(),std::numeric_limits<double>::max()))
399 ==ptToFill+(newNb+1)*nbOfComp)
400 std::copy(ptSrc+i*nbOfComp,ptSrc+(i+1)*nbOfComp,ptToFill+newNb*nbOfComp);
403 std::transform(ptSrc+i*nbOfComp,ptSrc+(i+1)*nbOfComp,ptToFill+newNb*nbOfComp,(double *)tmp,std::minus<double>());
404 std::transform((double *)tmp,((double *)tmp)+nbOfComp,(double *)tmp,std::ptr_fun<double,double>(fabs));
405 //if(!std::equal(ptSrc+i*nbOfComp,ptSrc+(i+1)*nbOfComp,ptToFill+newNb*nbOfComp))
406 if(*std::max_element((double *)tmp,((double *)tmp)+nbOfComp)>eps)
408 std::ostringstream oss;
409 oss << msg << " " << i << " and " << std::find(old2NewPtr,old2NewPtr+i,newNb)-old2NewPtr
410 << " have been merged and " << msg << " field on them are different !";
411 throw INTERP_KERNEL::Exception(oss.str().c_str());
418 void MEDCouplingFieldDiscretization::RenumberEntitiesFromN2OArr(const int *new2OldPtr, int new2OldSz, DataArrayDouble *arr, const char *msg)
420 int nbOfComp=arr->getNumberOfComponents();
421 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arrCpy=arr->deepCpy();
422 const double *ptSrc=arrCpy->getConstPointer();
423 arr->reAlloc(new2OldSz);
424 double *ptToFill=arr->getPointer();
425 for(int i=0;i<new2OldSz;i++)
427 int oldNb=new2OldPtr[i];
428 std::copy(ptSrc+oldNb*nbOfComp,ptSrc+(oldNb+1)*nbOfComp,ptToFill+i*nbOfComp);
432 MEDCouplingFieldDiscretization::~MEDCouplingFieldDiscretization()
436 TypeOfField MEDCouplingFieldDiscretizationP0::getEnum() const
442 * This method is simply called by MEDCouplingFieldDiscretization::deepCpy. It performs the deep copy of \a this.
444 * \sa MEDCouplingFieldDiscretization::deepCpy.
446 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationP0::clone() const
448 return new MEDCouplingFieldDiscretizationP0;
451 std::string MEDCouplingFieldDiscretizationP0::getStringRepr() const
453 return std::string(REPR);
456 const char *MEDCouplingFieldDiscretizationP0::getRepr() const
461 bool MEDCouplingFieldDiscretizationP0::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
465 reason="other spatial discretization is NULL, and this spatial discretization (P0) is defined.";
468 const MEDCouplingFieldDiscretizationP0 *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationP0 *>(other);
471 reason="Spatial discrtization of this is ON_CELLS, which is not the case of other.";
475 int MEDCouplingFieldDiscretizationP0::getNumberOfTuples(const MEDCouplingMesh *mesh) const throw(INTERP_KERNEL::Exception)
478 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getNumberOfTuples : NULL input mesh !");
479 return mesh->getNumberOfCells();
483 * mesh is not used here. It is not a bug !
485 int MEDCouplingFieldDiscretizationP0::getNumberOfTuplesExpectedRegardingCode(const MEDCouplingMesh *mesh, const std::vector<int>& code, const std::vector<const DataArrayInt *>& idsPerType) const throw(INTERP_KERNEL::Exception)
488 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getNumberOfTuplesExpectedRegardingCode : invalid input code !");
489 int nbOfSplit=(int)idsPerType.size();
490 int nbOfTypes=(int)code.size()/3;
492 for(int i=0;i<nbOfTypes;i++)
494 int nbOfEltInChunk=code[3*i+1];
496 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getNumberOfTuplesExpectedRegardingCode : invalid input code ! presence of negative value in a type !");
500 if(pos<0 || pos>=nbOfSplit)
502 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationP0::getNumberOfTuplesExpectedRegardingCode : input code points to pos " << pos << " in typeid " << i << " ! Should be in [0," << nbOfSplit << ") !";
503 throw INTERP_KERNEL::Exception(oss.str().c_str());
505 const DataArrayInt *ids(idsPerType[pos]);
506 if(!ids || !ids->isAllocated() || ids->getNumberOfComponents()!=1 || ids->getNumberOfTuples()!=nbOfEltInChunk || ids->getMinValueInArray()<0)
508 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationP0::getNumberOfTuplesExpectedRegardingCode : input pfl chunck at pos " << pos << " should have " << i << " tuples and one component and with ids all >=0 !";
509 throw INTERP_KERNEL::Exception(oss.str().c_str());
517 int MEDCouplingFieldDiscretizationP0::getNumberOfMeshPlaces(const MEDCouplingMesh *mesh) const
520 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getNumberOfMeshPlaces : NULL input mesh !");
521 return mesh->getNumberOfCells();
524 DataArrayInt *MEDCouplingFieldDiscretizationP0::getOffsetArr(const MEDCouplingMesh *mesh) const
527 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getOffsetArr : NULL input mesh !");
528 int nbOfTuples=mesh->getNumberOfCells();
529 DataArrayInt *ret=DataArrayInt::New();
530 ret->alloc(nbOfTuples+1,1);
535 void MEDCouplingFieldDiscretizationP0::renumberArraysForCell(const MEDCouplingMesh *mesh, const std::vector<DataArray *>& arrays,
536 const int *old2NewBg, bool check) throw(INTERP_KERNEL::Exception)
539 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::renumberArraysForCell : NULL input mesh !");
540 const int *array=old2NewBg;
542 array=DataArrayInt::CheckAndPreparePermutation(old2NewBg,old2NewBg+mesh->getNumberOfCells());
543 for(std::vector<DataArray *>::const_iterator it=arrays.begin();it!=arrays.end();it++)
546 (*it)->renumberInPlace(array);
549 free(const_cast<int *>(array));
552 DataArrayDouble *MEDCouplingFieldDiscretizationP0::getLocalizationOfDiscValues(const MEDCouplingMesh *mesh) const
555 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getLocalizationOfDiscValues : NULL input mesh !");
556 return mesh->getBarycenterAndOwner();
559 void MEDCouplingFieldDiscretizationP0::computeMeshRestrictionFromTupleIds(const MEDCouplingMesh *mesh, const int *tupleIdsBg, const int *tupleIdsEnd,
560 DataArrayInt *&cellRestriction, DataArrayInt *&trueTupleRestriction) const throw(INTERP_KERNEL::Exception)
563 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::computeMeshRestrictionFromTupleIds : NULL input mesh !");
564 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> tmp=DataArrayInt::New();
565 tmp->alloc((int)std::distance(tupleIdsBg,tupleIdsEnd),1);
566 std::copy(tupleIdsBg,tupleIdsEnd,tmp->getPointer());
567 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> tmp2(tmp->deepCpy());
568 cellRestriction=tmp.retn();
569 trueTupleRestriction=tmp2.retn();
572 void MEDCouplingFieldDiscretizationP0::reprQuickOverview(std::ostream& stream) const throw(INTERP_KERNEL::Exception)
574 stream << "P0 spatial discretization.";
577 void MEDCouplingFieldDiscretizationP0::checkCompatibilityWithNature(NatureOfField nat) const throw(INTERP_KERNEL::Exception)
581 void MEDCouplingFieldDiscretizationP0::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArray *da) const throw(INTERP_KERNEL::Exception)
584 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::checkCoherencyBetween : NULL input mesh or DataArray !");
585 if(mesh->getNumberOfCells()!=da->getNumberOfTuples())
587 std::ostringstream message;
588 message << "Field on cells invalid because there are " << mesh->getNumberOfCells();
589 message << " cells in mesh and " << da->getNumberOfTuples() << " tuples in field !";
590 throw INTERP_KERNEL::Exception(message.str().c_str());
594 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationP0::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
597 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getMeasureField : mesh instance specified is NULL !");
598 return mesh->getMeasureField(isAbs);
601 void MEDCouplingFieldDiscretizationP0::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
604 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getValueOn : NULL input mesh !");
605 int id=mesh->getCellContainingPoint(loc,_precision);
607 throw INTERP_KERNEL::Exception("Specified point is detected outside of mesh : unable to apply P0::getValueOn !");
608 arr->getTuple(id,res);
611 void MEDCouplingFieldDiscretizationP0::getValueOnPos(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, int i, int j, int k, double *res) const
613 const MEDCouplingCMesh *meshC=dynamic_cast<const MEDCouplingCMesh *>(mesh);
615 throw INTERP_KERNEL::Exception("P0::getValueOnPos is only accessible for structured meshes !");
616 int id=meshC->getCellIdFromPos(i,j,k);
617 arr->getTuple(id,res);
620 DataArrayDouble *MEDCouplingFieldDiscretizationP0::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const
623 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getValueOnMulti : NULL input mesh !");
624 std::vector<int> elts,eltsIndex;
625 mesh->getCellsContainingPoints(loc,nbOfPoints,_precision,elts,eltsIndex);
626 int spaceDim=mesh->getSpaceDimension();
627 int nbOfComponents=arr->getNumberOfComponents();
628 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
629 ret->alloc(nbOfPoints,nbOfComponents);
630 double *ptToFill=ret->getPointer();
631 for(int i=0;i<nbOfPoints;i++,ptToFill+=nbOfComponents)
632 if(eltsIndex[i+1]-eltsIndex[i]>=1)
633 arr->getTuple(elts[eltsIndex[i]],ptToFill);
636 std::ostringstream oss; oss << "Point #" << i << " with coordinates : (";
637 std::copy(loc+i*spaceDim,loc+(i+1)*spaceDim,std::ostream_iterator<double>(oss,", "));
638 oss << ") detected outside mesh : unable to apply P0::getValueOnMulti ! ";
639 throw INTERP_KERNEL::Exception(oss.str().c_str());
645 * Nothing to do. It's not a bug.
647 void MEDCouplingFieldDiscretizationP0::renumberValuesOnNodes(double , const int *, int newNbOfNodes, DataArrayDouble *) const
651 void MEDCouplingFieldDiscretizationP0::renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const
653 RenumberEntitiesFromO2NArr(epsOnVals,old2New,newSz,arr,"Cell");
656 void MEDCouplingFieldDiscretizationP0::renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const
658 RenumberEntitiesFromN2OArr(new2old,newSz,arr,"Cell");
662 * This method returns a tuple ids selection from cell ids selection [start;end).
663 * This method is called by MEDCouplingFieldDiscretizationP0::buildSubMeshData to return parameter \b di.
664 * Here for P0 it's very simple !
666 * \return a newly allocated array containing ids to select into the DataArrayDouble of the field.
669 DataArrayInt *MEDCouplingFieldDiscretizationP0::computeTupleIdsToSelectFromCellIds(const MEDCouplingMesh *mesh, const int *startCellIds, const int *endCellIds) const
671 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret=DataArrayInt::New();
672 ret->alloc((int)std::distance(startCellIds,endCellIds),1);
673 std::copy(startCellIds,endCellIds,ret->getPointer());
678 * This method returns a submesh of 'mesh' instance constituting cell ids contained in array defined as an interval [start;end).
679 * @param di is an array returned that specifies entity ids (here cells ids) in mesh 'mesh' of entity in returned submesh.
680 * Example : The first cell id of returned mesh has the (*di)[0] id in 'mesh'
682 * \sa MEDCouplingFieldDiscretizationP0::buildSubMeshDataRange
684 MEDCouplingMesh *MEDCouplingFieldDiscretizationP0::buildSubMeshData(const MEDCouplingMesh *mesh, const int *start, const int *end, DataArrayInt *&di) const
687 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::buildSubMeshData : NULL input mesh !");
688 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPart(start,end);
689 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> diSafe=DataArrayInt::New();
690 diSafe->alloc((int)std::distance(start,end),1);
691 std::copy(start,end,diSafe->getPointer());
697 * This method is strictly equivalent to MEDCouplingFieldDiscretizationP0::buildSubMeshData except that it is optimized for input defined as a range of cell ids.
699 * \param [out] beginOut Valid only if \a di is NULL
700 * \param [out] endOut Valid only if \a di is NULL
701 * \param [out] stepOut Valid only if \a di is NULL
702 * \param [out] di is an array returned that specifies entity ids (nodes, cells, Gauss points... ) in array if no output range is foundable.
704 * \sa MEDCouplingFieldDiscretizationP0::buildSubMeshData
706 MEDCouplingMesh *MEDCouplingFieldDiscretizationP0::buildSubMeshDataRange(const MEDCouplingMesh *mesh, int beginCellIds, int endCellIds, int stepCellIds, int& beginOut, int& endOut, int& stepOut, DataArrayInt *&di) const
709 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::buildSubMeshDataRange : NULL input mesh !");
710 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPartRange(beginCellIds,endCellIds,stepCellIds);
711 di=0; beginOut=beginCellIds; endOut=endCellIds; stepOut=stepCellIds;
715 int MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuples(const MEDCouplingMesh *mesh) const throw(INTERP_KERNEL::Exception)
718 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::getNumberOfTuples : NULL input mesh !");
719 return mesh->getNumberOfNodes();
723 * mesh is not used here. It is not a bug !
725 int MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuplesExpectedRegardingCode(const MEDCouplingMesh *mesh, const std::vector<int>& code, const std::vector<const DataArrayInt *>& idsPerType) const throw(INTERP_KERNEL::Exception)
728 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuplesExpectedRegardingCode : invalid input code !");
729 int nbOfSplit=(int)idsPerType.size();
730 int nbOfTypes=(int)code.size()/3;
732 for(int i=0;i<nbOfTypes;i++)
734 int nbOfEltInChunk=code[3*i+1];
736 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuplesExpectedRegardingCode : invalid input code ! presence of negative value in a type !");
740 if(pos<0 || pos>=nbOfSplit)
742 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuplesExpectedRegardingCode : input code points to pos " << pos << " in typeid " << i << " ! Should be in [0," << nbOfSplit << ") !";
743 throw INTERP_KERNEL::Exception(oss.str().c_str());
745 const DataArrayInt *ids(idsPerType[pos]);
746 if(!ids || !ids->isAllocated() || ids->getNumberOfComponents()!=1 || ids->getNumberOfTuples()!=nbOfEltInChunk || ids->getMinValueInArray()<0)
748 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuplesExpectedRegardingCode : input pfl chunck at pos " << pos << " should have " << i << " tuples and one component and with ids all >=0 !";
749 throw INTERP_KERNEL::Exception(oss.str().c_str());
757 int MEDCouplingFieldDiscretizationOnNodes::getNumberOfMeshPlaces(const MEDCouplingMesh *mesh) const
760 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::getNumberOfMeshPlaces : NULL input mesh !");
761 return mesh->getNumberOfNodes();
765 * Nothing to do here.
767 void MEDCouplingFieldDiscretizationOnNodes::renumberArraysForCell(const MEDCouplingMesh *, const std::vector<DataArray *>& arrays,
768 const int *old2NewBg, bool check) throw(INTERP_KERNEL::Exception)
772 DataArrayInt *MEDCouplingFieldDiscretizationOnNodes::getOffsetArr(const MEDCouplingMesh *mesh) const
775 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::getOffsetArr : NULL input mesh !");
776 int nbOfTuples=mesh->getNumberOfNodes();
777 DataArrayInt *ret=DataArrayInt::New();
778 ret->alloc(nbOfTuples+1,1);
783 DataArrayDouble *MEDCouplingFieldDiscretizationOnNodes::getLocalizationOfDiscValues(const MEDCouplingMesh *mesh) const
786 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::getLocalizationOfDiscValues : NULL input mesh !");
787 return mesh->getCoordinatesAndOwner();
790 void MEDCouplingFieldDiscretizationOnNodes::computeMeshRestrictionFromTupleIds(const MEDCouplingMesh *mesh, const int *tupleIdsBg, const int *tupleIdsEnd,
791 DataArrayInt *&cellRestriction, DataArrayInt *&trueTupleRestriction) const throw(INTERP_KERNEL::Exception)
794 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationOnNodes::computeMeshRestrictionFromTupleIds : NULL input mesh !");
795 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret1=mesh->getCellIdsFullyIncludedInNodeIds(tupleIdsBg,tupleIdsEnd);
796 const MEDCouplingUMesh *meshc=dynamic_cast<const MEDCouplingUMesh *>(mesh);
798 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationOnNodes::computeMeshRestrictionFromTupleIds : trying to subpart field on nodes by node ids ! Your mesh has to be unstructured !");
799 MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> meshPart=static_cast<MEDCouplingUMesh *>(meshc->buildPartOfMySelf(ret1->begin(),ret1->end(),true));
800 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret2=meshPart->computeFetchedNodeIds();
801 cellRestriction=ret1.retn();
802 trueTupleRestriction=ret2.retn();
805 void MEDCouplingFieldDiscretizationOnNodes::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArray *da) const throw(INTERP_KERNEL::Exception)
808 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::checkCoherencyBetween : NULL input mesh or DataArray !");
809 if(mesh->getNumberOfNodes()!=da->getNumberOfTuples())
811 std::ostringstream message;
812 message << "Field on nodes invalid because there are " << mesh->getNumberOfNodes();
813 message << " nodes in mesh and " << da->getNumberOfTuples() << " tuples in field !";
814 throw INTERP_KERNEL::Exception(message.str().c_str());
819 * This method returns a submesh of 'mesh' instance constituting cell ids contained in array defined as an interval [start;end).
820 * @param di is an array returned that specifies entity ids (here nodes ids) in mesh 'mesh' of entity in returned submesh.
821 * Example : The first node id of returned mesh has the (*di)[0] id in 'mesh'
823 MEDCouplingMesh *MEDCouplingFieldDiscretizationOnNodes::buildSubMeshData(const MEDCouplingMesh *mesh, const int *start, const int *end, DataArrayInt *&di) const
826 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::buildSubMeshData : NULL input mesh !");
827 DataArrayInt *diTmp=0;
828 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPartAndReduceNodes(start,end,diTmp);
829 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> diTmpSafe(diTmp);
830 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> di2=diTmpSafe->invertArrayO2N2N2O(ret->getNumberOfNodes());
836 * This method is strictly equivalent to MEDCouplingFieldDiscretizationNodes::buildSubMeshData except that it is optimized for input defined as a range of cell ids.
838 * \param [out] beginOut Valid only if \a di is NULL
839 * \param [out] endOut Valid only if \a di is NULL
840 * \param [out] stepOut Valid only if \a di is NULL
841 * \param [out] di is an array returned that specifies entity ids (nodes, cells, Gauss points... ) in array if no output range is foundable.
843 * \sa MEDCouplingFieldDiscretizationNodes::buildSubMeshData
845 MEDCouplingMesh *MEDCouplingFieldDiscretizationOnNodes::buildSubMeshDataRange(const MEDCouplingMesh *mesh, int beginCellIds, int endCellIds, int stepCellIds, int& beginOut, int& endOut, int& stepOut, DataArrayInt *&di) const
848 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationOnNodes::buildSubMeshDataRange : NULL input mesh !");
849 DataArrayInt *diTmp=0;
850 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPartRangeAndReduceNodes(beginCellIds,endCellIds,stepCellIds,beginOut,endOut,stepOut,diTmp);
853 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> diTmpSafe(diTmp);
854 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> di2=diTmpSafe->invertArrayO2N2N2O(ret->getNumberOfNodes());
861 * This method returns a tuple ids selection from cell ids selection [start;end).
862 * This method is called by MEDCouplingFieldDiscretizationOnNodes::buildSubMeshData to return parameter \b di.
863 * Here for P1 only nodes fetched by submesh of mesh[startCellIds:endCellIds) is returned !
865 * \return a newly allocated array containing ids to select into the DataArrayDouble of the field.
868 DataArrayInt *MEDCouplingFieldDiscretizationOnNodes::computeTupleIdsToSelectFromCellIds(const MEDCouplingMesh *mesh, const int *startCellIds, const int *endCellIds) const
871 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP1::computeTupleIdsToSelectFromCellIds : NULL input mesh !");
872 const MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> umesh=mesh->buildUnstructured();
873 MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> umesh2=static_cast<MEDCouplingUMesh *>(umesh->buildPartOfMySelf(startCellIds,endCellIds,true));
874 return umesh2->computeFetchedNodeIds();
877 void MEDCouplingFieldDiscretizationOnNodes::renumberValuesOnNodes(double epsOnVals, const int *old2NewPtr, int newNbOfNodes, DataArrayDouble *arr) const
879 RenumberEntitiesFromO2NArr(epsOnVals,old2NewPtr,newNbOfNodes,arr,"Node");
883 * Nothing to do it's not a bug.
885 void MEDCouplingFieldDiscretizationOnNodes::renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const
890 * Nothing to do it's not a bug.
892 void MEDCouplingFieldDiscretizationOnNodes::renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const
896 void MEDCouplingFieldDiscretizationOnNodes::getValueOnPos(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, int i, int j, int k, double *res) const
898 const MEDCouplingCMesh *meshC=dynamic_cast<const MEDCouplingCMesh *>(mesh);
900 throw INTERP_KERNEL::Exception("OnNodes::getValueOnPos(i,j,k) is only accessible for structured meshes !");
901 int id=meshC->getNodeIdFromPos(i,j,k);
902 arr->getTuple(id,res);
905 TypeOfField MEDCouplingFieldDiscretizationP1::getEnum() const
911 * This method is simply called by MEDCouplingFieldDiscretization::deepCpy. It performs the deep copy of \a this.
913 * \sa MEDCouplingFieldDiscretization::deepCpy.
915 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationP1::clone() const
917 return new MEDCouplingFieldDiscretizationP1;
920 std::string MEDCouplingFieldDiscretizationP1::getStringRepr() const
922 return std::string(REPR);
925 const char *MEDCouplingFieldDiscretizationP1::getRepr() const
930 bool MEDCouplingFieldDiscretizationP1::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
934 reason="other spatial discretization is NULL, and this spatial discretization (P1) is defined.";
937 const MEDCouplingFieldDiscretizationP1 *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationP1 *>(other);
940 reason="Spatial discrtization of this is ON_NODES, which is not the case of other.";
944 void MEDCouplingFieldDiscretizationP1::checkCompatibilityWithNature(NatureOfField nat) const throw(INTERP_KERNEL::Exception)
946 if(nat!=ConservativeVolumic)
947 throw INTERP_KERNEL::Exception("Invalid nature for P1 field : expected ConservativeVolumic !");
950 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationP1::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
953 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP1::getMeasureField : mesh instance specified is NULL !");
954 return mesh->getMeasureFieldOnNode(isAbs);
957 void MEDCouplingFieldDiscretizationP1::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
960 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP1::getValueOn : NULL input mesh !");
961 int id=mesh->getCellContainingPoint(loc,_precision);
963 throw INTERP_KERNEL::Exception("Specified point is detected outside of mesh : unable to apply P1::getValueOn !");
964 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(id);
965 if(type!=INTERP_KERNEL::NORM_SEG2 && type!=INTERP_KERNEL::NORM_TRI3 && type!=INTERP_KERNEL::NORM_TETRA4)
966 throw INTERP_KERNEL::Exception("P1 getValueOn is not specified for not simplex cells !");
967 getValueInCell(mesh,id,arr,loc,res);
971 * This method localizes a point defined by 'loc' in a cell with id 'cellId' into mesh 'mesh'.
972 * The result is put into res expected to be of size at least arr->getNumberOfComponents()
974 void MEDCouplingFieldDiscretizationP1::getValueInCell(const MEDCouplingMesh *mesh, int cellId, const DataArrayDouble *arr, const double *loc, double *res) const
977 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP1::getValueInCell : NULL input mesh !");
978 std::vector<int> conn;
979 std::vector<double> coo;
980 mesh->getNodeIdsOfCell(cellId,conn);
981 for(std::vector<int>::const_iterator iter=conn.begin();iter!=conn.end();iter++)
982 mesh->getCoordinatesOfNode(*iter,coo);
983 int spaceDim=mesh->getSpaceDimension();
984 std::size_t nbOfNodes=conn.size();
985 std::vector<const double *> vec(nbOfNodes);
986 for(std::size_t i=0;i<nbOfNodes;i++)
987 vec[i]=&coo[i*spaceDim];
988 INTERP_KERNEL::AutoPtr<double> tmp=new double[nbOfNodes];
989 INTERP_KERNEL::barycentric_coords(vec,loc,tmp);
990 int sz=arr->getNumberOfComponents();
991 INTERP_KERNEL::AutoPtr<double> tmp2=new double[sz];
992 std::fill(res,res+sz,0.);
993 for(std::size_t i=0;i<nbOfNodes;i++)
995 arr->getTuple(conn[i],(double *)tmp2);
996 std::transform((double *)tmp2,((double *)tmp2)+sz,(double *)tmp2,std::bind2nd(std::multiplies<double>(),tmp[i]));
997 std::transform(res,res+sz,(double *)tmp2,res,std::plus<double>());
1001 DataArrayDouble *MEDCouplingFieldDiscretizationP1::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const
1004 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP1::getValueOnMulti : NULL input mesh !");
1005 std::vector<int> elts,eltsIndex;
1006 mesh->getCellsContainingPoints(loc,nbOfPoints,_precision,elts,eltsIndex);
1007 int spaceDim=mesh->getSpaceDimension();
1008 int nbOfComponents=arr->getNumberOfComponents();
1009 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
1010 ret->alloc(nbOfPoints,nbOfComponents);
1011 double *ptToFill=ret->getPointer();
1012 for(int i=0;i<nbOfPoints;i++)
1013 if(eltsIndex[i+1]-eltsIndex[i]>=1)
1014 getValueInCell(mesh,elts[eltsIndex[i]],arr,loc+i*spaceDim,ptToFill+i*nbOfComponents);
1017 std::ostringstream oss; oss << "Point #" << i << " with coordinates : (";
1018 std::copy(loc+i*spaceDim,loc+(i+1)*spaceDim,std::ostream_iterator<double>(oss,", "));
1019 oss << ") detected outside mesh : unable to apply P1::getValueOnMulti ! ";
1020 throw INTERP_KERNEL::Exception(oss.str().c_str());
1025 void MEDCouplingFieldDiscretizationP1::reprQuickOverview(std::ostream& stream) const throw(INTERP_KERNEL::Exception)
1027 stream << "P1 spatial discretization.";
1030 MEDCouplingFieldDiscretizationPerCell::MEDCouplingFieldDiscretizationPerCell():_discr_per_cell(0)
1034 MEDCouplingFieldDiscretizationPerCell::~MEDCouplingFieldDiscretizationPerCell()
1037 _discr_per_cell->decrRef();
1041 * This constructor deep copies ParaMEDMEM::DataArrayInt instance from other (if any).
1043 MEDCouplingFieldDiscretizationPerCell::MEDCouplingFieldDiscretizationPerCell(const MEDCouplingFieldDiscretizationPerCell& other, const int *startCellIds, const int *endCellIds):_discr_per_cell(0)
1045 DataArrayInt *arr=other._discr_per_cell;
1048 if(startCellIds==0 && endCellIds==0)
1049 _discr_per_cell=arr->deepCpy();
1051 _discr_per_cell=arr->selectByTupleIdSafe(startCellIds,endCellIds);
1055 MEDCouplingFieldDiscretizationPerCell::MEDCouplingFieldDiscretizationPerCell(const MEDCouplingFieldDiscretizationPerCell& other, int beginCellIds, int endCellIds, int stepCellIds):_discr_per_cell(0)
1057 DataArrayInt *arr=other._discr_per_cell;
1060 _discr_per_cell=arr->selectByTupleId2(beginCellIds,endCellIds,stepCellIds);
1064 void MEDCouplingFieldDiscretizationPerCell::updateTime() const
1067 updateTimeWith(*_discr_per_cell);
1070 std::size_t MEDCouplingFieldDiscretizationPerCell::getHeapMemorySize() const
1074 ret+=_discr_per_cell->getHeapMemorySize();
1078 void MEDCouplingFieldDiscretizationPerCell::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArray *da) const throw(INTERP_KERNEL::Exception)
1080 if(!_discr_per_cell)
1081 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell has no discretization per cell !");
1083 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell::checkCoherencyBetween : NULL input mesh or DataArray !");
1084 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1085 if(nbOfTuples!=mesh->getNumberOfCells())
1086 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell has a discretization per cell but it's not matching the underlying mesh !");
1089 bool MEDCouplingFieldDiscretizationPerCell::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
1093 reason="other spatial discretization is NULL, and this spatial discretization (PerCell) is defined.";
1096 const MEDCouplingFieldDiscretizationPerCell *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationPerCell *>(other);
1099 reason="Spatial discretization of this is ON_GAUSS, which is not the case of other.";
1102 if(_discr_per_cell==0)
1103 return otherC->_discr_per_cell==0;
1104 if(otherC->_discr_per_cell==0)
1106 bool ret=_discr_per_cell->isEqualIfNotWhy(*otherC->_discr_per_cell,reason);
1108 reason.insert(0,"Field discretization per cell DataArrayInt given the discid per cell :");
1112 bool MEDCouplingFieldDiscretizationPerCell::isEqualWithoutConsideringStr(const MEDCouplingFieldDiscretization *other, double eps) const
1114 const MEDCouplingFieldDiscretizationPerCell *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationPerCell *>(other);
1117 if(_discr_per_cell==0)
1118 return otherC->_discr_per_cell==0;
1119 if(otherC->_discr_per_cell==0)
1121 return _discr_per_cell->isEqualWithoutConsideringStr(*otherC->_discr_per_cell);
1125 * This method is typically the first step of renumbering. The impact on _discr_per_cell is necessary here.
1126 * virtualy by this method.
1128 void MEDCouplingFieldDiscretizationPerCell::renumberCells(const int *old2NewBg, bool check) throw(INTERP_KERNEL::Exception)
1130 int nbCells=_discr_per_cell->getNumberOfTuples();
1131 const int *array=old2NewBg;
1133 array=DataArrayInt::CheckAndPreparePermutation(old2NewBg,old2NewBg+nbCells);
1135 DataArrayInt *dpc=_discr_per_cell->renumber(array);
1136 _discr_per_cell->decrRef();
1137 _discr_per_cell=dpc;
1140 free(const_cast<int *>(array));
1143 void MEDCouplingFieldDiscretizationPerCell::buildDiscrPerCellIfNecessary(const MEDCouplingMesh *mesh)
1146 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell::buildDiscrPerCellIfNecessary : NULL input mesh !");
1147 if(!_discr_per_cell)
1149 _discr_per_cell=DataArrayInt::New();
1150 int nbTuples=mesh->getNumberOfCells();
1151 _discr_per_cell->alloc(nbTuples,1);
1152 int *ptr=_discr_per_cell->getPointer();
1153 std::fill(ptr,ptr+nbTuples,DFT_INVALID_LOCID_VALUE);
1157 void MEDCouplingFieldDiscretizationPerCell::checkNoOrphanCells() const throw(INTERP_KERNEL::Exception)
1159 if(!_discr_per_cell)
1160 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell::checkNoOrphanCells : no discretization defined !");
1161 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> test=_discr_per_cell->getIdsEqual(DFT_INVALID_LOCID_VALUE);
1162 if(test->getNumberOfTuples()!=0)
1163 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell::checkNoOrphanCells : presence of orphan cells !");
1167 * This method is useful when 'this' describes a field discretization with several gauss discretization on a \b same cell type.
1168 * For example same NORM_TRI3 cells having 6 gauss points and others with 12 gauss points.
1169 * This method returns 2 arrays with same size : the return value and 'locIds' output parameter.
1170 * 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.
1171 * The return vector contains a set of newly created instance to deal with.
1172 * The returned vector represents a \b partition of cells ids with a gauss discretization set.
1174 * If no descretization is set in 'this' and exception will be thrown.
1176 std::vector<DataArrayInt *> MEDCouplingFieldDiscretizationPerCell::splitIntoSingleGaussDicrPerCellType(std::vector<int>& locIds) const throw(INTERP_KERNEL::Exception)
1178 if(!_discr_per_cell)
1179 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell::splitIntoSingleGaussDicrPerCellType : no descretization set !");
1180 return _discr_per_cell->partitionByDifferentValues(locIds);
1183 const DataArrayInt *MEDCouplingFieldDiscretizationPerCell::getArrayOfDiscIds() const
1185 return _discr_per_cell;
1188 void MEDCouplingFieldDiscretizationPerCell::setArrayOfDiscIds(const DataArrayInt *adids) throw(INTERP_KERNEL::Exception)
1190 if(adids!=_discr_per_cell)
1193 _discr_per_cell->decrRef();
1194 _discr_per_cell=const_cast<DataArrayInt *>(adids);
1196 _discr_per_cell->incrRef();
1201 MEDCouplingFieldDiscretizationGauss::MEDCouplingFieldDiscretizationGauss()
1205 MEDCouplingFieldDiscretizationGauss::MEDCouplingFieldDiscretizationGauss(const MEDCouplingFieldDiscretizationGauss& other, const int *startCellIds, const int *endCellIds):MEDCouplingFieldDiscretizationPerCell(other,startCellIds,endCellIds),_loc(other._loc)
1209 MEDCouplingFieldDiscretizationGauss::MEDCouplingFieldDiscretizationGauss(const MEDCouplingFieldDiscretizationGauss& other, int beginCellIds, int endCellIds, int stepCellIds):MEDCouplingFieldDiscretizationPerCell(other,beginCellIds,endCellIds,stepCellIds),_loc(other._loc)
1213 TypeOfField MEDCouplingFieldDiscretizationGauss::getEnum() const
1218 bool MEDCouplingFieldDiscretizationGauss::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
1222 reason="other spatial discretization is NULL, and this spatial discretization (Gauss) is defined.";
1225 const MEDCouplingFieldDiscretizationGauss *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationGauss *>(other);
1228 reason="Spatial discrtization of this is ON_GAUSS, which is not the case of other.";
1231 if(!MEDCouplingFieldDiscretizationPerCell::isEqualIfNotWhy(other,eps,reason))
1233 if(_loc.size()!=otherC->_loc.size())
1235 reason="Gauss spatial discretization : localization sizes differ";
1238 std::size_t sz=_loc.size();
1239 for(std::size_t i=0;i<sz;i++)
1240 if(!_loc[i].isEqual(otherC->_loc[i],eps))
1242 std::ostringstream oss; oss << "Gauss spatial discretization : Localization #" << i << " differ from this to other.";
1249 bool MEDCouplingFieldDiscretizationGauss::isEqualWithoutConsideringStr(const MEDCouplingFieldDiscretization *other, double eps) const
1251 const MEDCouplingFieldDiscretizationGauss *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationGauss *>(other);
1254 if(!MEDCouplingFieldDiscretizationPerCell::isEqualWithoutConsideringStr(other,eps))
1256 if(_loc.size()!=otherC->_loc.size())
1258 std::size_t sz=_loc.size();
1259 for(std::size_t i=0;i<sz;i++)
1260 if(!_loc[i].isEqual(otherC->_loc[i],eps))
1266 * This method is simply called by MEDCouplingFieldDiscretization::deepCpy. It performs the deep copy of \a this.
1268 * \sa MEDCouplingFieldDiscretization::deepCpy.
1270 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationGauss::clone() const
1272 return new MEDCouplingFieldDiscretizationGauss(*this);
1275 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationGauss::clonePart(const int *startCellIds, const int *endCellIds) const
1277 return new MEDCouplingFieldDiscretizationGauss(*this,startCellIds,endCellIds);
1280 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationGauss::clonePartRange(int beginCellIds, int endCellIds, int stepCellIds) const
1282 return new MEDCouplingFieldDiscretizationGauss(*this,beginCellIds,endCellIds,stepCellIds);
1285 std::string MEDCouplingFieldDiscretizationGauss::getStringRepr() const
1287 std::ostringstream oss; oss << REPR << "." << std::endl;
1290 if(_discr_per_cell->isAllocated())
1292 oss << "Discretization per cell : ";
1293 std::copy(_discr_per_cell->begin(),_discr_per_cell->end(),std::ostream_iterator<int>(oss,", "));
1297 oss << "Presence of " << _loc.size() << " localizations." << std::endl;
1299 for(std::vector<MEDCouplingGaussLocalization>::const_iterator it=_loc.begin();it!=_loc.end();it++,i++)
1301 oss << "+++++ Localization #" << i << " +++++" << std::endl;
1302 oss << (*it).getStringRepr();
1303 oss << "++++++++++" << std::endl;
1308 std::size_t MEDCouplingFieldDiscretizationGauss::getHeapMemorySize() const
1310 std::size_t ret=_loc.capacity()*sizeof(MEDCouplingGaussLocalization);
1311 for(std::vector<MEDCouplingGaussLocalization>::const_iterator it=_loc.begin();it!=_loc.end();it++)
1312 ret+=(*it).getHeapMemorySize();
1313 return MEDCouplingFieldDiscretizationPerCell::getHeapMemorySize()+ret;
1316 const char *MEDCouplingFieldDiscretizationGauss::getRepr() const
1322 * mesh is not used here. It is not a bug !
1324 int MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode(const MEDCouplingMesh *mesh, const std::vector<int>& code, const std::vector<const DataArrayInt *>& idsPerType) const throw(INTERP_KERNEL::Exception)
1326 if(!_discr_per_cell || !_discr_per_cell->isAllocated() || _discr_per_cell->getNumberOfComponents()!=1)
1327 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode");
1328 if(code.size()%3!=0)
1329 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode : invalid input code !");
1330 int nbOfSplit=(int)idsPerType.size();
1331 int nbOfTypes=(int)code.size()/3;
1333 for(int i=0;i<nbOfTypes;i++)
1335 int nbOfEltInChunk=code[3*i+1];
1336 if(nbOfEltInChunk<0)
1337 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode : invalid input code ! presence of negative value in a type !");
1338 int pos=code[3*i+2];
1341 if(pos<0 || pos>=nbOfSplit)
1343 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode : input code points to pos " << pos << " in typeid " << i << " ! Should be in [0," << nbOfSplit << ") !";
1344 throw INTERP_KERNEL::Exception(oss.str().c_str());
1346 const DataArrayInt *ids(idsPerType[pos]);
1347 if(!ids || !ids->isAllocated() || ids->getNumberOfComponents()!=1 || ids->getNumberOfTuples()!=nbOfEltInChunk || ids->getMinValueInArray()<0)
1349 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode : input pfl chunck at pos " << pos << " should have " << i << " tuples and one component and with ids all >=0 !";
1350 throw INTERP_KERNEL::Exception(oss.str().c_str());
1353 ret+=nbOfEltInChunk;
1355 if(ret!=_discr_per_cell->getNumberOfTuples())
1357 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode : input code points to " << ret << " cells whereas discretization percell array lgth is " << _discr_per_cell->getNumberOfTuples() << " !";
1359 return getNumberOfTuples(0);//0 is not an error ! It is to be sure that input mesh is not used
1362 int MEDCouplingFieldDiscretizationGauss::getNumberOfTuples(const MEDCouplingMesh *) const throw(INTERP_KERNEL::Exception)
1365 if (_discr_per_cell == 0)
1366 throw INTERP_KERNEL::Exception("Discretization is not initialized!");
1367 const int *dcPtr=_discr_per_cell->getConstPointer();
1368 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1369 int maxSz=(int)_loc.size();
1370 for(const int *w=dcPtr;w!=dcPtr+nbOfTuples;w++)
1372 if(*w>=0 && *w<maxSz)
1373 ret+=_loc[*w].getNumberOfGaussPt();
1376 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getNumberOfTuples : At cell #" << std::distance(dcPtr,w) << " localization id is " << *w << " should be in [0," << maxSz << ") !";
1377 throw INTERP_KERNEL::Exception(oss.str().c_str());
1383 int MEDCouplingFieldDiscretizationGauss::getNumberOfMeshPlaces(const MEDCouplingMesh *mesh) const
1386 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getNumberOfMeshPlaces : NULL input mesh !");
1387 return mesh->getNumberOfCells();
1391 * This method is redevelopped for performance reasons, but it is equivalent to a call to MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField
1392 * and a call to DataArrayDouble::computeOffsets2 on the returned array.
1394 DataArrayInt *MEDCouplingFieldDiscretizationGauss::getOffsetArr(const MEDCouplingMesh *mesh) const
1397 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getOffsetArr : NULL input mesh !");
1398 int nbOfTuples=mesh->getNumberOfCells();
1399 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret=DataArrayInt::New();
1400 ret->alloc(nbOfTuples+1,1);
1401 int *retPtr=ret->getPointer();
1402 const int *start=_discr_per_cell->getConstPointer();
1403 if(_discr_per_cell->getNumberOfTuples()!=nbOfTuples)
1404 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getOffsetArr : mismatch between the mesh and the discretization ids array length !");
1405 int maxPossible=(int)_loc.size();
1407 for(int i=0;i<nbOfTuples;i++,start++)
1409 if(*start>=0 && *start<maxPossible)
1410 retPtr[i+1]=retPtr[i]+_loc[*start].getNumberOfGaussPt();
1413 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getOffsetArr : At position #" << i << " the locid = " << *start << " whereas it should be in [0," << maxPossible << ") !";
1414 throw INTERP_KERNEL::Exception(oss.str().c_str());
1420 void MEDCouplingFieldDiscretizationGauss::renumberArraysForCell(const MEDCouplingMesh *mesh, const std::vector<DataArray *>& arrays,
1421 const int *old2NewBg, bool check) throw(INTERP_KERNEL::Exception)
1424 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::renumberArraysForCell : NULL input mesh !");
1425 const int *array=old2NewBg;
1427 array=DataArrayInt::CheckAndPreparePermutation(old2NewBg,old2NewBg+mesh->getNumberOfCells());
1428 int nbOfCells=_discr_per_cell->getNumberOfTuples();
1429 int nbOfTuples=getNumberOfTuples(0);
1430 const int *dcPtr=_discr_per_cell->getConstPointer();
1431 int *array2=new int[nbOfTuples];//stores the final conversion array old2New to give to arrays in renumberInPlace.
1432 int *array3=new int[nbOfCells];//store for each cell in present dcp array (already renumbered) the offset needed by each cell in new numbering.
1434 for(int i=1;i<nbOfCells;i++)
1435 array3[i]=array3[i-1]+_loc[dcPtr[i-1]].getNumberOfGaussPt();
1437 for(int i=0;i<nbOfCells;i++)
1439 int nbOfGaussPt=_loc[dcPtr[array[i]]].getNumberOfGaussPt();
1440 for(int k=0;k<nbOfGaussPt;k++,j++)
1441 array2[j]=array3[array[i]]+k;
1444 for(std::vector<DataArray *>::const_iterator it=arrays.begin();it!=arrays.end();it++)
1446 (*it)->renumberInPlace(array2);
1449 free(const_cast<int*>(array));
1452 DataArrayDouble *MEDCouplingFieldDiscretizationGauss::getLocalizationOfDiscValues(const MEDCouplingMesh *mesh) const
1455 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getLocalizationOfDiscValues : NULL input mesh !");
1456 checkNoOrphanCells();
1457 MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> umesh=mesh->buildUnstructured();//in general do nothing
1458 int nbOfTuples=getNumberOfTuples(mesh);
1459 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
1460 int spaceDim=mesh->getSpaceDimension();
1461 ret->alloc(nbOfTuples,spaceDim);
1462 std::vector< int > locIds;
1463 std::vector<DataArrayInt *> parts=splitIntoSingleGaussDicrPerCellType(locIds);
1464 std::vector< MEDCouplingAutoRefCountObjectPtr<DataArrayInt> > parts2(parts.size());
1465 std::copy(parts.begin(),parts.end(),parts2.begin());
1466 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> offsets=buildNbOfGaussPointPerCellField();
1467 offsets->computeOffsets();
1468 const int *ptrOffsets=offsets->getConstPointer();
1469 const double *coords=umesh->getCoords()->getConstPointer();
1470 const int *connI=umesh->getNodalConnectivityIndex()->getConstPointer();
1471 const int *conn=umesh->getNodalConnectivity()->getConstPointer();
1472 double *valsToFill=ret->getPointer();
1473 for(std::size_t i=0;i<parts2.size();i++)
1475 INTERP_KERNEL::GaussCoords calculator;
1477 const MEDCouplingGaussLocalization& cli=_loc[locIds[i]];//curLocInfo
1478 INTERP_KERNEL::NormalizedCellType typ=cli.getType();
1479 const std::vector<double>& wg=cli.getWeights();
1480 calculator.addGaussInfo(typ,INTERP_KERNEL::CellModel::GetCellModel(typ).getDimension(),
1481 &cli.getGaussCoords()[0],(int)wg.size(),&cli.getRefCoords()[0],
1482 INTERP_KERNEL::CellModel::GetCellModel(typ).getNumberOfNodes());
1484 int nbt=parts2[i]->getNumberOfTuples();
1485 for(const int *w=parts2[i]->getConstPointer();w!=parts2[i]->getConstPointer()+nbt;w++)
1486 calculator.calculateCoords(cli.getType(),coords,spaceDim,conn+connI[*w]+1,valsToFill+spaceDim*(ptrOffsets[*w]));
1488 ret->copyStringInfoFrom(*umesh->getCoords());
1492 void MEDCouplingFieldDiscretizationGauss::computeMeshRestrictionFromTupleIds(const MEDCouplingMesh *mesh, const int *tupleIdsBg, const int *tupleIdsEnd,
1493 DataArrayInt *&cellRestriction, DataArrayInt *&trueTupleRestriction) const throw(INTERP_KERNEL::Exception)
1496 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::computeMeshRestrictionFromTupleIds : NULL input mesh !");
1497 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> tmp=DataArrayInt::New(); tmp->alloc((int)std::distance(tupleIdsBg,tupleIdsEnd),1);
1498 std::copy(tupleIdsBg,tupleIdsEnd,tmp->getPointer());
1500 tmp=tmp->buildUnique();
1501 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=buildNbOfGaussPointPerCellField();
1502 nbOfNodesPerCell->computeOffsets2();
1503 nbOfNodesPerCell->searchRangesInListOfIds(tmp,cellRestriction,trueTupleRestriction);
1509 void MEDCouplingFieldDiscretizationGauss::checkCompatibilityWithNature(NatureOfField nat) const throw(INTERP_KERNEL::Exception)
1513 void MEDCouplingFieldDiscretizationGauss::getTinySerializationIntInformation(std::vector<int>& tinyInfo) const
1517 val=_discr_per_cell->getNumberOfTuples();
1518 tinyInfo.push_back(val);
1519 tinyInfo.push_back((int)_loc.size());
1521 tinyInfo.push_back(-1);
1523 tinyInfo.push_back(_loc[0].getDimension());
1524 for(std::vector<MEDCouplingGaussLocalization>::const_iterator iter=_loc.begin();iter!=_loc.end();iter++)
1525 (*iter).pushTinySerializationIntInfo(tinyInfo);
1528 void MEDCouplingFieldDiscretizationGauss::getTinySerializationDbleInformation(std::vector<double>& tinyInfo) const
1530 for(std::vector<MEDCouplingGaussLocalization>::const_iterator iter=_loc.begin();iter!=_loc.end();iter++)
1531 (*iter).pushTinySerializationDblInfo(tinyInfo);
1534 void MEDCouplingFieldDiscretizationGauss::getSerializationIntArray(DataArrayInt *& arr) const
1538 arr=_discr_per_cell;
1541 void MEDCouplingFieldDiscretizationGauss::resizeForUnserialization(const std::vector<int>& tinyInfo, DataArrayInt *& arr)
1543 int val=tinyInfo[0];
1546 _discr_per_cell=DataArrayInt::New();
1547 _discr_per_cell->alloc(val,1);
1551 arr=_discr_per_cell;
1552 int nbOfLoc=tinyInfo[1];
1554 int dim=tinyInfo[2];
1557 delta=((int)tinyInfo.size()-3)/nbOfLoc;
1558 for(int i=0;i<nbOfLoc;i++)
1560 std::vector<int> tmp(tinyInfo.begin()+3+i*delta,tinyInfo.begin()+3+(i+1)*delta);
1561 MEDCouplingGaussLocalization elt=MEDCouplingGaussLocalization::BuildNewInstanceFromTinyInfo(dim,tmp);
1562 _loc.push_back(elt);
1566 void MEDCouplingFieldDiscretizationGauss::finishUnserialization(const std::vector<double>& tinyInfo)
1568 double *tmp=new double[tinyInfo.size()];
1569 std::copy(tinyInfo.begin(),tinyInfo.end(),tmp);
1570 const double *work=tmp;
1571 for(std::vector<MEDCouplingGaussLocalization>::iterator iter=_loc.begin();iter!=_loc.end();iter++)
1572 work=(*iter).fillWithValues(work);
1576 double MEDCouplingFieldDiscretizationGauss::getIJK(const MEDCouplingMesh *mesh, const DataArrayDouble *da,
1577 int cellId, int nodeIdInCell, int compoId) const throw(INTERP_KERNEL::Exception)
1579 int offset=getOffsetOfCell(cellId);
1580 return da->getIJ(offset+nodeIdInCell,compoId);
1583 void MEDCouplingFieldDiscretizationGauss::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArray *da) const throw(INTERP_KERNEL::Exception)
1586 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::checkCoherencyBetween : NULL input mesh or DataArray !");
1587 MEDCouplingFieldDiscretizationPerCell::checkCoherencyBetween(mesh,da);
1588 for(std::vector<MEDCouplingGaussLocalization>::const_iterator iter=_loc.begin();iter!=_loc.end();iter++)
1589 (*iter).checkCoherency();
1590 int nbOfDesc=(int)_loc.size();
1591 int nbOfCells=mesh->getNumberOfCells();
1592 const int *dc=_discr_per_cell->getConstPointer();
1593 for(int i=0;i<nbOfCells;i++)
1597 std::ostringstream oss; oss << "Cell # " << i << " of mesh \"" << mesh->getName() << "\" has an undefined gauss location ! Should never happend !";
1598 throw INTERP_KERNEL::Exception(oss.str().c_str());
1602 std::ostringstream oss; oss << "Cell # " << i << " of mesh \"" << mesh->getName() << "\" has no gauss location !";
1603 throw INTERP_KERNEL::Exception(oss.str().c_str());
1605 if(mesh->getTypeOfCell(i)!=_loc[dc[i]].getType())
1607 std::ostringstream oss; oss << "Types of mesh and gauss location mismatch for cell # " << i;
1608 throw INTERP_KERNEL::Exception(oss.str().c_str());
1611 int nbOfTuples=getNumberOfTuples(mesh);
1612 if(nbOfTuples!=da->getNumberOfTuples())
1614 std::ostringstream oss; oss << "Invalid number of tuples in the array : expecting " << nbOfTuples << " !";
1615 throw INTERP_KERNEL::Exception(oss.str().c_str());
1619 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationGauss::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
1622 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getMeasureField : mesh instance specified is NULL !");
1623 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=mesh->getMeasureField(isAbs);
1624 const double *volPtr=vol->getArray()->begin();
1625 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> ret=MEDCouplingFieldDouble::New(ON_GAUSS_PT);
1627 ret->setDiscretization(const_cast<MEDCouplingFieldDiscretizationGauss *>(this));
1628 if(!_discr_per_cell)
1629 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getMeasureField : no discr per cell array not defined ! spatial localization is incorrect !");
1630 _discr_per_cell->checkAllocated();
1631 if(_discr_per_cell->getNumberOfComponents()!=1)
1632 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getMeasureField : no discr per cell array defined but with nb of components different from 1 !");
1633 if(_discr_per_cell->getNumberOfTuples()!=vol->getNumberOfTuples())
1634 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getMeasureField : no discr per cell array defined but mismatch between nb of cells of mesh and size of spatial disr array !");
1635 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> offset=getOffsetArr(mesh);
1636 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arr=DataArrayDouble::New(); arr->alloc(getNumberOfTuples(mesh),1);
1638 double *arrPtr=arr->getPointer();
1639 const int *offsetPtr=offset->getConstPointer();
1640 int maxGaussLoc=(int)_loc.size();
1641 std::vector<int> locIds;
1642 std::vector<DataArrayInt *> ids=splitIntoSingleGaussDicrPerCellType(locIds);
1643 std::vector< MEDCouplingAutoRefCountObjectPtr<DataArrayInt> > ids2(ids.size()); std::copy(ids.begin(),ids.end(),ids2.begin());
1644 for(std::size_t i=0;i<locIds.size();i++)
1646 const DataArrayInt *curIds=ids[i];
1647 int locId=locIds[i];
1648 if(locId>=0 && locId<maxGaussLoc)
1650 const MEDCouplingGaussLocalization& loc=_loc[locId];
1651 int nbOfGaussPt=loc.getNumberOfGaussPt();
1652 INTERP_KERNEL::AutoPtr<double> weights=new double[nbOfGaussPt];
1653 double sum=std::accumulate(loc.getWeights().begin(),loc.getWeights().end(),0.);
1654 std::transform(loc.getWeights().begin(),loc.getWeights().end(),(double *)weights,std::bind2nd(std::multiplies<double>(),1./sum));
1655 for(const int *cellId=curIds->begin();cellId!=curIds->end();cellId++)
1656 for(int j=0;j<nbOfGaussPt;j++)
1657 arrPtr[offsetPtr[*cellId]+j]=weights[j]*volPtr[*cellId];
1661 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getMeasureField : Presence of localization id " << locId << " in cell #" << curIds->getIJ(0,0) << " ! Must be in [0," << maxGaussLoc << ") !";
1662 throw INTERP_KERNEL::Exception(oss.str().c_str());
1665 ret->synchronizeTimeWithSupport();
1669 void MEDCouplingFieldDiscretizationGauss::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
1671 throw INTERP_KERNEL::Exception("Not implemented yet !");
1674 void MEDCouplingFieldDiscretizationGauss::getValueOnPos(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, int i, int j, int k, double *res) const
1676 throw INTERP_KERNEL::Exception("getValueOnPos(i,j,k) : Not applyable for Gauss points !");
1679 DataArrayDouble *MEDCouplingFieldDiscretizationGauss::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const
1681 throw INTERP_KERNEL::Exception("getValueOnMulti : Not implemented yet for gauss points !");
1684 MEDCouplingMesh *MEDCouplingFieldDiscretizationGauss::buildSubMeshData(const MEDCouplingMesh *mesh, const int *start, const int *end, DataArrayInt *&di) const
1687 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::buildSubMeshData : NULL input mesh !");
1688 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> diSafe=computeTupleIdsToSelectFromCellIds(mesh,start,end);
1689 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPart(start,end);
1695 * This method is strictly equivalent to MEDCouplingFieldDiscretizationGauss::buildSubMeshData except that it is optimized for input defined as a range of cell ids.
1697 * \param [out] beginOut Valid only if \a di is NULL
1698 * \param [out] endOut Valid only if \a di is NULL
1699 * \param [out] stepOut Valid only if \a di is NULL
1700 * \param [out] di is an array returned that specifies entity ids (nodes, cells, Gauss points... ) in array if no output range is foundable.
1702 * \sa MEDCouplingFieldDiscretizationGauss::buildSubMeshData
1704 MEDCouplingMesh *MEDCouplingFieldDiscretizationGauss::buildSubMeshDataRange(const MEDCouplingMesh *mesh, int beginCellIds, int endCellIds, int stepCellIds, int& beginOut, int& endOut, int& stepOut, DataArrayInt *&di) const
1706 if(stepCellIds!=1)//even for stepCellIds==-1 the output will not be a range
1707 return MEDCouplingFieldDiscretization::buildSubMeshDataRange(mesh,beginCellIds,endCellIds,stepCellIds,beginOut,endOut,stepOut,di);
1709 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::buildSubMeshDataRange : NULL input mesh !");
1710 if(!_discr_per_cell)
1711 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::buildSubMeshDataRange : no discretization array set !");
1712 di=0; beginOut=0; endOut=0; stepOut=stepCellIds;
1713 const char msg[]="MEDCouplingFieldDiscretizationGauss::buildSubMeshDataRange : cell #";
1714 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1715 const int *w=_discr_per_cell->begin();
1716 int nbMaxOfLocId=(int)_loc.size();
1717 for(int i=0;i<nbOfTuples;i++,w++)
1719 if(*w!=DFT_INVALID_LOCID_VALUE)
1721 if(*w>=0 && *w<nbMaxOfLocId)
1723 int delta=_loc[*w].getNumberOfGaussPt();
1731 { std::ostringstream oss; oss << msg << i << " has invalid id (" << *w << ") ! Should be in [0," << nbMaxOfLocId << ") !"; throw INTERP_KERNEL::Exception(oss.str().c_str()); }
1734 { std::ostringstream oss; oss << msg << i << " is detected as orphan !"; throw INTERP_KERNEL::Exception(oss.str().c_str()); }
1736 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPartRange(beginCellIds,endCellIds,stepCellIds);
1741 * This method returns a tuple ids selection from cell ids selection [start;end).
1742 * This method is called by MEDCouplingFieldDiscretizationGauss::buildSubMeshData to return parameter \b di.
1744 * \return a newly allocated array containing ids to select into the DataArrayDouble of the field.
1747 DataArrayInt *MEDCouplingFieldDiscretizationGauss::computeTupleIdsToSelectFromCellIds(const MEDCouplingMesh *mesh, const int *startCellIds, const int *endCellIds) const
1750 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::computeTupleIdsToSelectFromCellIds : null mesh !");
1751 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=buildNbOfGaussPointPerCellField();//check of _discr_per_cell not NULL pointer
1752 int nbOfCells=mesh->getNumberOfCells();
1753 if(_discr_per_cell->getNumberOfTuples()!=nbOfCells)
1754 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::computeTupleIdsToSelectFromCellIds : mismatch of nb of tuples of cell ids array and number of cells !");
1755 nbOfNodesPerCell->computeOffsets2();
1756 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> sel=DataArrayInt::New(); sel->useArray(startCellIds,false,CPP_DEALLOC,(int)std::distance(startCellIds,endCellIds),1);
1757 return sel->buildExplicitArrByRanges(nbOfNodesPerCell);
1761 * No implementation needed !
1763 void MEDCouplingFieldDiscretizationGauss::renumberValuesOnNodes(double , const int *, int newNbOfNodes, DataArrayDouble *) const
1767 void MEDCouplingFieldDiscretizationGauss::renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const
1769 throw INTERP_KERNEL::Exception("Not implemented yet !");
1772 void MEDCouplingFieldDiscretizationGauss::renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const
1774 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 !");
1777 void MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnType(const MEDCouplingMesh *mesh, INTERP_KERNEL::NormalizedCellType type, const std::vector<double>& refCoo,
1778 const std::vector<double>& gsCoo, const std::vector<double>& wg) throw(INTERP_KERNEL::Exception)
1781 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnType : NULL input mesh !");
1782 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
1783 if((int)cm.getDimension()!=mesh->getMeshDimension())
1785 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnType : mismatch of dimensions ! MeshDim==" << mesh->getMeshDimension();
1786 oss << " whereas Type '" << cm.getRepr() << "' has dimension " << cm.getDimension() << " !";
1787 throw INTERP_KERNEL::Exception(oss.str().c_str());
1789 buildDiscrPerCellIfNecessary(mesh);
1790 int id=(int)_loc.size();
1791 MEDCouplingGaussLocalization elt(type,refCoo,gsCoo,wg);
1792 _loc.push_back(elt);
1793 int *ptr=_discr_per_cell->getPointer();
1794 int nbCells=mesh->getNumberOfCells();
1795 for(int i=0;i<nbCells;i++)
1796 if(mesh->getTypeOfCell(i)==type)
1798 zipGaussLocalizations();
1801 void MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnCells(const MEDCouplingMesh *mesh, const int *begin, const int *end, const std::vector<double>& refCoo,
1802 const std::vector<double>& gsCoo, const std::vector<double>& wg) throw(INTERP_KERNEL::Exception)
1805 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnCells : NULL input mesh !");
1806 buildDiscrPerCellIfNecessary(mesh);
1807 if(std::distance(begin,end)<1)
1808 throw INTERP_KERNEL::Exception("Size of [begin,end) must be equal or greater than 1 !");
1809 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(*begin);
1810 MEDCouplingGaussLocalization elt(type,refCoo,gsCoo,wg);
1811 int id=(int)_loc.size();
1812 int *ptr=_discr_per_cell->getPointer();
1813 for(const int *w=begin+1;w!=end;w++)
1815 if(mesh->getTypeOfCell(*w)!=type)
1817 std::ostringstream oss; oss << "The cell with id " << *w << " has been detected to be incompatible in the [begin,end) array specified !";
1818 throw INTERP_KERNEL::Exception(oss.str().c_str());
1822 for(const int *w2=begin;w2!=end;w2++)
1825 _loc.push_back(elt);
1826 zipGaussLocalizations();
1829 void MEDCouplingFieldDiscretizationGauss::clearGaussLocalizations() throw(INTERP_KERNEL::Exception)
1833 _discr_per_cell->decrRef();
1839 void MEDCouplingFieldDiscretizationGauss::setGaussLocalization(int locId, const MEDCouplingGaussLocalization& loc) throw(INTERP_KERNEL::Exception)
1842 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::setGaussLocalization : localization id has to be >=0 !");
1843 int sz=(int)_loc.size();
1844 MEDCouplingGaussLocalization gLoc(INTERP_KERNEL::NORM_ERROR);
1846 _loc.resize(locId+1,gLoc);
1850 void MEDCouplingFieldDiscretizationGauss::resizeLocalizationVector(int newSz) throw(INTERP_KERNEL::Exception)
1853 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::resizeLocalizationVector : new size has to be >=0 !");
1854 MEDCouplingGaussLocalization gLoc(INTERP_KERNEL::NORM_ERROR);
1855 _loc.resize(newSz,gLoc);
1858 MEDCouplingGaussLocalization& MEDCouplingFieldDiscretizationGauss::getGaussLocalization(int locId) throw(INTERP_KERNEL::Exception)
1860 checkLocalizationId(locId);
1864 int MEDCouplingFieldDiscretizationGauss::getNbOfGaussLocalization() const throw(INTERP_KERNEL::Exception)
1866 return (int)_loc.size();
1869 int MEDCouplingFieldDiscretizationGauss::getGaussLocalizationIdOfOneCell(int cellId) const throw(INTERP_KERNEL::Exception)
1871 if(!_discr_per_cell)
1872 throw INTERP_KERNEL::Exception("No Gauss localization still set !");
1873 int locId=_discr_per_cell->begin()[cellId];
1875 throw INTERP_KERNEL::Exception("No Gauss localization set for the specified cell !");
1879 int MEDCouplingFieldDiscretizationGauss::getGaussLocalizationIdOfOneType(INTERP_KERNEL::NormalizedCellType type) const throw(INTERP_KERNEL::Exception)
1881 std::set<int> ret=getGaussLocalizationIdsOfOneType(type);
1883 throw INTERP_KERNEL::Exception("No gauss discretization found for the specified type !");
1885 throw INTERP_KERNEL::Exception("Several gauss discretizations have been found for the specified type !");
1886 return *ret.begin();
1889 std::set<int> MEDCouplingFieldDiscretizationGauss::getGaussLocalizationIdsOfOneType(INTERP_KERNEL::NormalizedCellType type) const throw(INTERP_KERNEL::Exception)
1891 if(!_discr_per_cell)
1892 throw INTERP_KERNEL::Exception("No Gauss localization still set !");
1895 for(std::vector<MEDCouplingGaussLocalization>::const_iterator iter=_loc.begin();iter!=_loc.end();iter++,id++)
1896 if((*iter).getType()==type)
1901 void MEDCouplingFieldDiscretizationGauss::getCellIdsHavingGaussLocalization(int locId, std::vector<int>& cellIds) const throw(INTERP_KERNEL::Exception)
1903 if(locId<0 || locId>=(int)_loc.size())
1904 throw INTERP_KERNEL::Exception("Invalid locId given : must be in range [0:getNbOfGaussLocalization()) !");
1905 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1906 const int *ptr=_discr_per_cell->getConstPointer();
1907 for(int i=0;i<nbOfTuples;i++)
1909 cellIds.push_back(i);
1912 const MEDCouplingGaussLocalization& MEDCouplingFieldDiscretizationGauss::getGaussLocalization(int locId) const throw(INTERP_KERNEL::Exception)
1914 checkLocalizationId(locId);
1918 void MEDCouplingFieldDiscretizationGauss::checkLocalizationId(int locId) const throw(INTERP_KERNEL::Exception)
1920 if(locId<0 || locId>=(int)_loc.size())
1921 throw INTERP_KERNEL::Exception("Invalid locId given : must be in range [0:getNbOfGaussLocalization()) !");
1924 int MEDCouplingFieldDiscretizationGauss::getOffsetOfCell(int cellId) const throw(INTERP_KERNEL::Exception)
1927 const int *start=_discr_per_cell->getConstPointer();
1928 for(const int *w=start;w!=start+cellId;w++)
1929 ret+=_loc[*w].getNumberOfGaussPt();
1934 * This method do the assumption that there is no orphan cell. If there is an exception is thrown.
1935 * This method makes the assumption too that '_discr_per_cell' is defined. If not an exception is thrown.
1936 * This method returns a newly created array with number of tuples equals to '_discr_per_cell->getNumberOfTuples' and number of components equal to 1.
1937 * The i_th tuple in returned array is the number of gauss point if the corresponding cell.
1939 DataArrayInt *MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField() const throw(INTERP_KERNEL::Exception)
1941 if(!_discr_per_cell)
1942 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField : no discretization array set !");
1943 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1944 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret=DataArrayInt::New();
1945 const int *w=_discr_per_cell->begin();
1946 ret->alloc(nbOfTuples,1);
1947 int *valsToFill=ret->getPointer();
1948 int nbMaxOfLocId=(int)_loc.size();
1949 for(int i=0;i<nbOfTuples;i++,w++)
1950 if(*w!=DFT_INVALID_LOCID_VALUE)
1952 if(*w>=0 && *w<nbMaxOfLocId)
1953 valsToFill[i]=_loc[*w].getNumberOfGaussPt();
1956 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField : cell #" << i << " has invalid id (" << *w << ") ! Should be in [0," << nbMaxOfLocId << ") !";
1957 throw INTERP_KERNEL::Exception(oss.str().c_str());
1962 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField : cell #" << i << " is detected as orphan !";
1963 throw INTERP_KERNEL::Exception(oss.str().c_str());
1968 void MEDCouplingFieldDiscretizationGauss::reprQuickOverview(std::ostream& stream) const throw(INTERP_KERNEL::Exception)
1970 stream << "Gauss points spatial discretization.";
1974 * This method makes the assumption that _discr_per_cell is set.
1975 * This method reduces as much as possible number size of _loc.
1976 * This method is useful when several set on same cells has been done and that some Gauss Localization are no more used.
1978 void MEDCouplingFieldDiscretizationGauss::zipGaussLocalizations()
1980 const int *start=_discr_per_cell->begin();
1981 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1982 INTERP_KERNEL::AutoPtr<int> tmp=new int[_loc.size()];
1983 std::fill((int *)tmp,(int *)tmp+_loc.size(),-2);
1984 for(const int *w=start;w!=start+nbOfTuples;w++)
1988 for(int i=0;i<(int)_loc.size();i++)
1991 if(fid==(int)_loc.size())
1994 int *start2=_discr_per_cell->getPointer();
1995 for(int *w2=start2;w2!=start2+nbOfTuples;w2++)
1998 std::vector<MEDCouplingGaussLocalization> tmpLoc;
1999 for(int i=0;i<(int)_loc.size();i++)
2001 tmpLoc.push_back(_loc[tmp[i]]);
2005 MEDCouplingFieldDiscretizationGaussNE::MEDCouplingFieldDiscretizationGaussNE()
2009 TypeOfField MEDCouplingFieldDiscretizationGaussNE::getEnum() const
2015 * This method is simply called by MEDCouplingFieldDiscretization::deepCpy. It performs the deep copy of \a this.
2017 * \sa MEDCouplingFieldDiscretization::deepCpy.
2019 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationGaussNE::clone() const
2021 return new MEDCouplingFieldDiscretizationGaussNE(*this);
2024 std::string MEDCouplingFieldDiscretizationGaussNE::getStringRepr() const
2026 return std::string(REPR);
2029 const char *MEDCouplingFieldDiscretizationGaussNE::getRepr() const
2034 bool MEDCouplingFieldDiscretizationGaussNE::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
2038 reason="other spatial discretization is NULL, and this spatial discretization (GaussNE) is defined.";
2041 const MEDCouplingFieldDiscretizationGaussNE *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationGaussNE *>(other);
2044 reason="Spatial discrtization of this is ON_GAUSS_NE, which is not the case of other.";
2048 int MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode(const MEDCouplingMesh *mesh, const std::vector<int>& code, const std::vector<const DataArrayInt *>& idsPerType) const throw(INTERP_KERNEL::Exception)
2050 if(code.size()%3!=0)
2051 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode : invalid input code !");
2052 int nbOfSplit=(int)idsPerType.size();
2053 int nbOfTypes=(int)code.size()/3;
2055 for(int i=0;i<nbOfTypes;i++)
2057 int nbOfEltInChunk=code[3*i+1];
2058 if(nbOfEltInChunk<0)
2059 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode : invalid input code ! presence of negative value in a type !");
2060 int pos=code[3*i+2];
2063 if(pos<0 || pos>=nbOfSplit)
2065 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode : input code points to pos " << pos << " in typeid " << i << " ! Should be in [0," << nbOfSplit << ") !";
2066 throw INTERP_KERNEL::Exception(oss.str().c_str());
2068 const DataArrayInt *ids(idsPerType[pos]);
2069 if(!ids || !ids->isAllocated() || ids->getNumberOfComponents()!=1 || ids->getNumberOfTuples()!=nbOfEltInChunk || ids->getMinValueInArray()<0)
2071 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode : input pfl chunck at pos " << pos << " should have " << i << " tuples and one component and with ids all >=0 !";
2072 throw INTERP_KERNEL::Exception(oss.str().c_str());
2075 ret+=nbOfEltInChunk;
2078 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode : NULL input mesh !");
2079 if(ret!=mesh->getNumberOfCells())
2081 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode : input code points to " << ret << " number of cells should be " << mesh->getNumberOfCells() << " !";
2083 return getNumberOfTuples(mesh);
2086 int MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuples(const MEDCouplingMesh *mesh) const throw(INTERP_KERNEL::Exception)
2089 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuples : NULL input mesh !");
2091 int nbOfCells=mesh->getNumberOfCells();
2092 for(int i=0;i<nbOfCells;i++)
2094 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
2095 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2097 throw INTERP_KERNEL::Exception("Not implemented yet Gauss node on elements for polygons and polyedrons !");
2098 ret+=cm.getNumberOfNodes();
2103 int MEDCouplingFieldDiscretizationGaussNE::getNumberOfMeshPlaces(const MEDCouplingMesh *mesh) const
2106 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getNumberOfMeshPlaces : NULL input mesh !");
2107 return mesh->getNumberOfCells();
2110 DataArrayInt *MEDCouplingFieldDiscretizationGaussNE::getOffsetArr(const MEDCouplingMesh *mesh) const
2113 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getOffsetArr : NULL input mesh !");
2114 int nbOfTuples=mesh->getNumberOfCells();
2115 DataArrayInt *ret=DataArrayInt::New();
2116 ret->alloc(nbOfTuples+1,1);
2117 int *retPtr=ret->getPointer();
2119 for(int i=0;i<nbOfTuples;i++)
2121 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
2122 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2124 throw INTERP_KERNEL::Exception("Not implemented yet Gauss node on elements for polygons and polyedrons !");
2125 retPtr[i+1]=retPtr[i]+cm.getNumberOfNodes();
2130 void MEDCouplingFieldDiscretizationGaussNE::renumberArraysForCell(const MEDCouplingMesh *mesh, const std::vector<DataArray *>& arrays,
2131 const int *old2NewBg, bool check) throw(INTERP_KERNEL::Exception)
2134 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::renumberArraysForCell : NULL input mesh !");
2135 const int *array=old2NewBg;
2137 array=DataArrayInt::CheckAndPreparePermutation(old2NewBg,old2NewBg+mesh->getNumberOfCells());
2138 int nbOfCells=mesh->getNumberOfCells();
2139 int nbOfTuples=getNumberOfTuples(mesh);
2140 int *array2=new int[nbOfTuples];//stores the final conversion array old2New to give to arrays in renumberInPlace.
2141 int *array3=new int[nbOfCells];//store for each cell in after renumbering the offset needed by each cell in new numbering.
2143 for(int i=1;i<nbOfCells;i++)
2145 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell((int)std::distance(array,std::find(array,array+nbOfCells,i-1)));
2146 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2147 array3[i]=array3[i-1]+cm.getNumberOfNodes();
2150 for(int i=0;i<nbOfCells;i++)
2152 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
2153 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2154 for(int k=0;k<(int)cm.getNumberOfNodes();k++,j++)
2155 array2[j]=array3[array[i]]+k;
2158 for(std::vector<DataArray *>::const_iterator it=arrays.begin();it!=arrays.end();it++)
2160 (*it)->renumberInPlace(array2);
2163 free(const_cast<int *>(array));
2166 DataArrayDouble *MEDCouplingFieldDiscretizationGaussNE::getLocalizationOfDiscValues(const MEDCouplingMesh *mesh) const
2169 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getLocalizationOfDiscValues : NULL input mesh !");
2170 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
2171 MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> umesh=mesh->buildUnstructured();//in general do nothing
2172 int nbOfTuples=getNumberOfTuples(umesh);
2173 int spaceDim=mesh->getSpaceDimension();
2174 ret->alloc(nbOfTuples,spaceDim);
2175 const double *coords=umesh->getCoords()->begin();
2176 const int *connI=umesh->getNodalConnectivityIndex()->getConstPointer();
2177 const int *conn=umesh->getNodalConnectivity()->getConstPointer();
2178 int nbCells=umesh->getNumberOfCells();
2179 double *retPtr=ret->getPointer();
2180 for(int i=0;i<nbCells;i++,connI++)
2181 for(const int *w=conn+connI[0]+1;w!=conn+connI[1];w++)
2183 retPtr=std::copy(coords+(*w)*spaceDim,coords+((*w)+1)*spaceDim,retPtr);
2188 * Reimplemented from MEDCouplingFieldDiscretization::integral for performance reason. The default implementation is valid too for GAUSS_NE spatial discretization.
2190 void MEDCouplingFieldDiscretizationGaussNE::integral(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, bool isWAbs, double *res) const throw(INTERP_KERNEL::Exception)
2193 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::integral : input mesh or array is null !");
2194 int nbOfCompo=arr->getNumberOfComponents();
2195 std::fill(res,res+nbOfCompo,0.);
2197 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=mesh->getMeasureField(isWAbs);
2198 std::set<INTERP_KERNEL::NormalizedCellType> types=mesh->getAllGeoTypes();
2199 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=mesh->computeNbOfNodesPerCell();
2200 nbOfNodesPerCell->computeOffsets2();
2201 const double *arrPtr=arr->begin(),*volPtr=vol->getArray()->begin();
2202 for(std::set<INTERP_KERNEL::NormalizedCellType>::const_iterator it=types.begin();it!=types.end();it++)
2204 std::size_t wArrSz=-1;
2205 const double *wArr=GetWeightArrayFromGeometricType(*it,wArrSz);
2206 INTERP_KERNEL::AutoPtr<double> wArr2=new double[wArrSz];
2207 double sum=std::accumulate(wArr,wArr+wArrSz,0.);
2208 std::transform(wArr,wArr+wArrSz,(double *)wArr2,std::bind2nd(std::multiplies<double>(),1./sum));
2209 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ids=mesh->giveCellsWithType(*it);
2210 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ids2=ids->buildExplicitArrByRanges(nbOfNodesPerCell);
2211 const int *ptIds2=ids2->begin(),*ptIds=ids->begin();
2212 int nbOfCellsWithCurGeoType=ids->getNumberOfTuples();
2213 for(int i=0;i<nbOfCellsWithCurGeoType;i++,ptIds++,ptIds2+=wArrSz)
2215 for(int k=0;k<nbOfCompo;k++)
2218 for(std::size_t j=0;j<wArrSz;j++)
2219 tmp+=arrPtr[nbOfCompo*ptIds2[j]+k]*wArr2[j];
2220 res[k]+=tmp*volPtr[*ptIds];
2226 const double *MEDCouplingFieldDiscretizationGaussNE::GetWeightArrayFromGeometricType(INTERP_KERNEL::NormalizedCellType geoType, std::size_t& lgth) throw(INTERP_KERNEL::Exception)
2230 case INTERP_KERNEL::NORM_SEG2:
2231 lgth=(int)sizeof(FGP_SEG2)/sizeof(double);
2233 case INTERP_KERNEL::NORM_SEG3:
2234 lgth=(int)sizeof(FGP_SEG3)/sizeof(double);
2236 case INTERP_KERNEL::NORM_SEG4:
2237 lgth=(int)sizeof(FGP_SEG4)/sizeof(double);
2239 case INTERP_KERNEL::NORM_TRI3:
2240 lgth=(int)sizeof(FGP_TRI3)/sizeof(double);
2242 case INTERP_KERNEL::NORM_TRI6:
2243 lgth=(int)sizeof(FGP_TRI6)/sizeof(double);
2245 case INTERP_KERNEL::NORM_TRI7:
2246 lgth=(int)sizeof(FGP_TRI7)/sizeof(double);
2248 case INTERP_KERNEL::NORM_QUAD4:
2249 lgth=(int)sizeof(FGP_QUAD4)/sizeof(double);
2251 case INTERP_KERNEL::NORM_QUAD9:
2252 lgth=(int)sizeof(FGP_QUAD9)/sizeof(double);
2254 case INTERP_KERNEL::NORM_TETRA4:
2255 lgth=(int)sizeof(FGP_TETRA4)/sizeof(double);
2257 case INTERP_KERNEL::NORM_PENTA6:
2258 lgth=(int)sizeof(FGP_PENTA6)/sizeof(double);
2260 case INTERP_KERNEL::NORM_HEXA8:
2261 lgth=(int)sizeof(FGP_HEXA8)/sizeof(double);
2263 case INTERP_KERNEL::NORM_HEXA27:
2264 lgth=(int)sizeof(FGP_HEXA27)/sizeof(double);
2266 case INTERP_KERNEL::NORM_PYRA5:
2267 lgth=(int)sizeof(FGP_PYRA5)/sizeof(double);
2270 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::GetWeightArrayFromGeometricType : only SEG[2,3,4], TRI[3,6,7], QUAD[4,9], TETRA4, PENTA6, HEXA[8,27], PYRA5 supported !");
2274 const double *MEDCouplingFieldDiscretizationGaussNE::GetRefCoordsFromGeometricType(INTERP_KERNEL::NormalizedCellType geoType, std::size_t& lgth) throw(INTERP_KERNEL::Exception)
2278 case INTERP_KERNEL::NORM_SEG2:
2279 lgth=(int)sizeof(REF_SEG2)/sizeof(double);
2281 case INTERP_KERNEL::NORM_SEG3:
2282 lgth=(int)sizeof(REF_SEG3)/sizeof(double);
2284 case INTERP_KERNEL::NORM_SEG4:
2285 lgth=(int)sizeof(REF_SEG4)/sizeof(double);
2287 case INTERP_KERNEL::NORM_TRI3:
2288 lgth=(int)sizeof(REF_TRI3)/sizeof(double);
2290 case INTERP_KERNEL::NORM_TRI6:
2291 lgth=(int)sizeof(REF_TRI6)/sizeof(double);
2293 case INTERP_KERNEL::NORM_TRI7:
2294 lgth=(int)sizeof(REF_TRI7)/sizeof(double);
2296 case INTERP_KERNEL::NORM_QUAD4:
2297 lgth=(int)sizeof(REF_QUAD4)/sizeof(double);
2299 case INTERP_KERNEL::NORM_QUAD8:
2300 lgth=(int)sizeof(REF_QUAD8)/sizeof(double);
2302 case INTERP_KERNEL::NORM_QUAD9:
2303 lgth=(int)sizeof(REF_QUAD9)/sizeof(double);
2305 case INTERP_KERNEL::NORM_TETRA4:
2306 lgth=(int)sizeof(REF_TETRA4)/sizeof(double);
2308 case INTERP_KERNEL::NORM_TETRA10:
2309 lgth=(int)sizeof(REF_TETRA10)/sizeof(double);
2311 case INTERP_KERNEL::NORM_PENTA6:
2312 lgth=(int)sizeof(REF_PENTA6)/sizeof(double);
2314 case INTERP_KERNEL::NORM_PENTA15:
2315 lgth=(int)sizeof(REF_PENTA15)/sizeof(double);
2317 case INTERP_KERNEL::NORM_HEXA8:
2318 lgth=(int)sizeof(REF_HEXA8)/sizeof(double);
2320 case INTERP_KERNEL::NORM_HEXA20:
2321 lgth=(int)sizeof(REF_HEXA20)/sizeof(double);
2323 case INTERP_KERNEL::NORM_HEXA27:
2324 lgth=(int)sizeof(REF_HEXA27)/sizeof(double);
2326 case INTERP_KERNEL::NORM_PYRA5:
2327 lgth=(int)sizeof(REF_PYRA5)/sizeof(double);
2329 case INTERP_KERNEL::NORM_PYRA13:
2330 lgth=(int)sizeof(REF_PYRA13)/sizeof(double);
2333 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::GetRefCoordsFromGeometricType : only SEG[2,3,4], TRI[3,6,7], QUAD[4,8,9], TETRA[4,10], PENTA[6,15], HEXA[8,20,27], PYRA[5,13] supported !");
2337 void MEDCouplingFieldDiscretizationGaussNE::computeMeshRestrictionFromTupleIds(const MEDCouplingMesh *mesh, const int *tupleIdsBg, const int *tupleIdsEnd,
2338 DataArrayInt *&cellRestriction, DataArrayInt *&trueTupleRestriction) const throw(INTERP_KERNEL::Exception)
2341 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::computeMeshRestrictionFromTupleIds : NULL input mesh !");
2342 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> tmp=DataArrayInt::New(); tmp->alloc((int)std::distance(tupleIdsBg,tupleIdsEnd),1);
2343 std::copy(tupleIdsBg,tupleIdsEnd,tmp->getPointer());
2345 tmp=tmp->buildUnique();
2346 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=mesh->computeNbOfNodesPerCell();
2347 nbOfNodesPerCell->computeOffsets2();
2348 nbOfNodesPerCell->searchRangesInListOfIds(tmp,cellRestriction,trueTupleRestriction);
2351 void MEDCouplingFieldDiscretizationGaussNE::checkCompatibilityWithNature(NatureOfField nat) const throw(INTERP_KERNEL::Exception)
2355 double MEDCouplingFieldDiscretizationGaussNE::getIJK(const MEDCouplingMesh *mesh, const DataArrayDouble *da,
2356 int cellId, int nodeIdInCell, int compoId) const throw(INTERP_KERNEL::Exception)
2359 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getIJK : NULL input mesh !");
2361 for(int i=0;i<cellId;i++)
2363 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
2364 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2365 offset+=cm.getNumberOfNodes();
2367 return da->getIJ(offset+nodeIdInCell,compoId);
2370 void MEDCouplingFieldDiscretizationGaussNE::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArray *da) const throw(INTERP_KERNEL::Exception)
2372 int nbOfTuples=getNumberOfTuples(mesh);
2373 if(nbOfTuples!=da->getNumberOfTuples())
2375 std::ostringstream oss; oss << "Invalid number of tuples in the array : expecting " << nbOfTuples << " !";
2376 throw INTERP_KERNEL::Exception(oss.str().c_str());
2380 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationGaussNE::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
2383 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getMeasureField : mesh instance specified is NULL !");
2384 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=mesh->getMeasureField(isAbs);
2385 const double *volPtr=vol->getArray()->begin();
2386 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> ret=MEDCouplingFieldDouble::New(ON_GAUSS_NE);
2389 std::set<INTERP_KERNEL::NormalizedCellType> types=mesh->getAllGeoTypes();
2390 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=mesh->computeNbOfNodesPerCell();
2391 int nbTuples=nbOfNodesPerCell->accumulate(0);
2392 nbOfNodesPerCell->computeOffsets2();
2393 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arr=DataArrayDouble::New(); arr->alloc(nbTuples,1);
2395 double *arrPtr=arr->getPointer();
2396 for(std::set<INTERP_KERNEL::NormalizedCellType>::const_iterator it=types.begin();it!=types.end();it++)
2398 std::size_t wArrSz=-1;
2399 const double *wArr=GetWeightArrayFromGeometricType(*it,wArrSz);
2400 INTERP_KERNEL::AutoPtr<double> wArr2=new double[wArrSz];
2401 double sum=std::accumulate(wArr,wArr+wArrSz,0.);
2402 std::transform(wArr,wArr+wArrSz,(double *)wArr2,std::bind2nd(std::multiplies<double>(),1./sum));
2403 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ids=mesh->giveCellsWithType(*it);
2404 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ids2=ids->buildExplicitArrByRanges(nbOfNodesPerCell);
2405 const int *ptIds2=ids2->begin(),*ptIds=ids->begin();
2406 int nbOfCellsWithCurGeoType=ids->getNumberOfTuples();
2407 for(int i=0;i<nbOfCellsWithCurGeoType;i++,ptIds++)
2408 for(std::size_t j=0;j<wArrSz;j++,ptIds2++)
2409 arrPtr[*ptIds2]=wArr2[j]*volPtr[*ptIds];
2411 ret->synchronizeTimeWithSupport();
2415 void MEDCouplingFieldDiscretizationGaussNE::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
2417 throw INTERP_KERNEL::Exception("Not implemented yet !");
2420 void MEDCouplingFieldDiscretizationGaussNE::getValueOnPos(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, int i, int j, int k, double *res) const
2422 throw INTERP_KERNEL::Exception("getValueOnPos(i,j,k) : Not applyable for Gauss points !");
2425 DataArrayDouble *MEDCouplingFieldDiscretizationGaussNE::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const
2427 throw INTERP_KERNEL::Exception("getValueOnMulti : Not implemented for Gauss NE !");
2430 MEDCouplingMesh *MEDCouplingFieldDiscretizationGaussNE::buildSubMeshData(const MEDCouplingMesh *mesh, const int *start, const int *end, DataArrayInt *&di) const
2433 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::buildSubMeshData : NULL input mesh !");
2434 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> diSafe=computeTupleIdsToSelectFromCellIds(mesh,start,end);
2435 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPart(start,end);
2441 * This method is strictly equivalent to MEDCouplingFieldDiscretizationGauss::buildSubMeshData except that it is optimized for input defined as a range of cell ids.
2443 * \param [out] beginOut Valid only if \a di is NULL
2444 * \param [out] endOut Valid only if \a di is NULL
2445 * \param [out] stepOut Valid only if \a di is NULL
2446 * \param [out] di is an array returned that specifies entity ids (nodes, cells, Gauss points... ) in array if no output range is foundable.
2448 * \sa MEDCouplingFieldDiscretizationGauss::buildSubMeshData
2450 MEDCouplingMesh *MEDCouplingFieldDiscretizationGaussNE::buildSubMeshDataRange(const MEDCouplingMesh *mesh, int beginCellIds, int endCellIds, int stepCellIds, int& beginOut, int& endOut, int& stepOut, DataArrayInt *&di) const
2452 if(stepCellIds!=1)//even for stepCellIds==-1 the output will not be a range
2453 return MEDCouplingFieldDiscretization::buildSubMeshDataRange(mesh,beginCellIds,endCellIds,stepCellIds,beginOut,endOut,stepOut,di);
2455 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::buildSubMeshDataRange : NULL input mesh !");
2456 int nbOfCells=mesh->getNumberOfCells();
2457 di=0; beginOut=0; endOut=0; stepOut=stepCellIds;
2458 const char msg[]="MEDCouplingFieldDiscretizationGaussNE::buildSubMeshDataRange : cell #";
2459 for(int i=0;i<nbOfCells;i++)
2461 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
2462 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2464 { std::ostringstream oss; oss << msg << i << " presence of dynamic cell (polygons and polyedrons) ! Not implemented !"; throw INTERP_KERNEL::Exception(oss.str().c_str()); }
2465 int delta=cm.getNumberOfNodes();
2472 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPartRange(beginCellIds,endCellIds,stepCellIds);
2478 * This method returns a tuple ids selection from cell ids selection [start;end).
2479 * This method is called by MEDCouplingFieldDiscretizationGaussNE::buildSubMeshData to return parameter \b di.
2481 * \return a newly allocated array containing ids to select into the DataArrayDouble of the field.
2484 DataArrayInt *MEDCouplingFieldDiscretizationGaussNE::computeTupleIdsToSelectFromCellIds(const MEDCouplingMesh *mesh, const int *startCellIds, const int *endCellIds) const
2487 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::computeTupleIdsToSelectFromCellIds : null mesh !");
2488 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=mesh->computeNbOfNodesPerCell();
2489 nbOfNodesPerCell->computeOffsets2();
2490 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> sel=DataArrayInt::New(); sel->useArray(startCellIds,false,CPP_DEALLOC,(int)std::distance(startCellIds,endCellIds),1);
2491 return sel->buildExplicitArrByRanges(nbOfNodesPerCell);
2495 * No implementation needed !
2497 void MEDCouplingFieldDiscretizationGaussNE::renumberValuesOnNodes(double , const int *, int newNbOfNodes, DataArrayDouble *) const
2501 void MEDCouplingFieldDiscretizationGaussNE::renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const
2503 throw INTERP_KERNEL::Exception("Not implemented yet !");
2506 void MEDCouplingFieldDiscretizationGaussNE::renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const
2508 throw INTERP_KERNEL::Exception("Not implemented yet !");
2511 void MEDCouplingFieldDiscretizationGaussNE::reprQuickOverview(std::ostream& stream) const throw(INTERP_KERNEL::Exception)
2513 stream << "Gauss points on nodes per element spatial discretization.";
2516 MEDCouplingFieldDiscretizationGaussNE::MEDCouplingFieldDiscretizationGaussNE(const MEDCouplingFieldDiscretizationGaussNE& other):MEDCouplingFieldDiscretization(other)
2520 TypeOfField MEDCouplingFieldDiscretizationKriging::getEnum() const
2525 const char *MEDCouplingFieldDiscretizationKriging::getRepr() const
2531 * This method is simply called by MEDCouplingFieldDiscretization::deepCpy. It performs the deep copy of \a this.
2533 * \sa MEDCouplingFieldDiscretization::deepCpy.
2535 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationKriging::clone() const
2537 return new MEDCouplingFieldDiscretizationKriging;
2540 std::string MEDCouplingFieldDiscretizationKriging::getStringRepr() const
2542 return std::string(REPR);
2545 void MEDCouplingFieldDiscretizationKriging::checkCompatibilityWithNature(NatureOfField nat) const throw(INTERP_KERNEL::Exception)
2547 if(nat!=ConservativeVolumic)
2548 throw INTERP_KERNEL::Exception("Invalid nature for Kriging field : expected ConservativeVolumic !");
2551 bool MEDCouplingFieldDiscretizationKriging::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
2555 reason="other spatial discretization is NULL, and this spatial discretization (Kriginig) is defined.";
2558 const MEDCouplingFieldDiscretizationKriging *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationKriging *>(other);
2561 reason="Spatial discrtization of this is ON_NODES_KR, which is not the case of other.";
2565 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationKriging::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
2568 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationKriging::getMeasureField : mesh instance specified is NULL !");
2569 throw INTERP_KERNEL::Exception("getMeasureField on FieldDiscretizationKriging : not implemented yet !");
2572 void MEDCouplingFieldDiscretizationKriging::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
2574 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> res2=MEDCouplingFieldDiscretizationKriging::getValueOnMulti(arr,mesh,loc,1);
2575 std::copy(res2->begin(),res2->end(),res);
2578 DataArrayDouble *MEDCouplingFieldDiscretizationKriging::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfTargetPoints) const
2581 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationKriging::getValueOnMulti : NULL input mesh !");
2582 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coords=getLocalizationOfDiscValues(mesh);
2583 int nbOfPts=coords->getNumberOfTuples();
2584 int dimension=coords->getNumberOfComponents();
2587 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> KnewiK=computeVectorOfCoefficients(mesh,arr,delta);
2589 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> locArr=DataArrayDouble::New();
2590 locArr->useArray(loc,false,CPP_DEALLOC,nbOfTargetPoints,dimension);
2591 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrix2=coords->buildEuclidianDistanceDenseMatrixWith(locArr);
2592 operateOnDenseMatrix(mesh->getSpaceDimension(),nbOfPts*nbOfTargetPoints,matrix2->getPointer());
2593 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrix3=DataArrayDouble::New();
2594 matrix3->alloc((nbOfPts+delta)*nbOfTargetPoints,1);
2595 double *work=matrix3->getPointer();
2596 const double *workCst=matrix2->getConstPointer();
2597 const double *workCst2=loc;
2598 for(int i=0;i<nbOfTargetPoints;i++,workCst+=nbOfPts,workCst2+=delta-1)
2600 for(int j=0;j<nbOfPts;j++)
2601 work[j*nbOfTargetPoints+i]=workCst[j];
2602 work[nbOfPts*nbOfTargetPoints+i]=1.0;
2603 for(int j=0;j<delta-1;j++)
2604 work[(nbOfPts+1+j)*nbOfTargetPoints+i]=workCst2[j];
2607 int nbOfCompo=arr->getNumberOfComponents();
2608 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
2609 ret->alloc(nbOfTargetPoints,nbOfCompo);
2610 INTERP_KERNEL::matrixProduct(KnewiK->getConstPointer(),1,nbOfPts+delta,matrix3->getConstPointer(),nbOfPts+delta,nbOfTargetPoints*nbOfCompo,ret->getPointer());
2614 void MEDCouplingFieldDiscretizationKriging::reprQuickOverview(std::ostream& stream) const throw(INTERP_KERNEL::Exception)
2616 stream << "Kriging spatial discretization.";
2620 * 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
2621 * number of tuples should be equal to the number of representing points in \a mesh.
2623 * \param [in] mesh is the sources of nodes on which kriging will be done regarding the parameters and the value of \c this->getSpaceDimension()
2624 * \param [in] arr input field DataArrayDouble whose number of tuples must be equal to the number of nodes in \a mesh
2625 * \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.
2626 * Whatever the value of \a isDrift the number of tuples of returned DataArrayDouble will be equal to \c arr->getNumberOfTuples() + \a isDrift.
2627 * \return a newly allocated array containing coefficients including or not drift coefficient at the end depending the value of \a isDrift parameter.
2629 DataArrayDouble *MEDCouplingFieldDiscretizationKriging::computeVectorOfCoefficients(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, int& isDrift) const
2632 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationKriging::computeVectorOfCoefficients : NULL input mesh !");
2633 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coords=getLocalizationOfDiscValues(mesh);
2634 int nbOfPts=coords->getNumberOfTuples();
2635 //int dimension=coords->getNumberOfComponents();
2636 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrix=coords->buildEuclidianDistanceDenseMatrix();
2637 operateOnDenseMatrix(mesh->getSpaceDimension(),nbOfPts*nbOfPts,matrix->getPointer());
2639 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrixWithDrift=performDrift(matrix,coords,isDrift);
2640 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrixInv=DataArrayDouble::New();
2641 matrixInv->alloc((nbOfPts+isDrift)*(nbOfPts+isDrift),1);
2642 INTERP_KERNEL::inverseMatrix(matrixWithDrift->getConstPointer(),nbOfPts+isDrift,matrixInv->getPointer());
2644 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> KnewiK=DataArrayDouble::New();
2645 KnewiK->alloc((nbOfPts+isDrift)*1,1);
2646 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arr2=DataArrayDouble::New();
2647 arr2->alloc((nbOfPts+isDrift)*1,1);
2648 double *work=std::copy(arr->begin(),arr->end(),arr2->getPointer());
2649 std::fill(work,work+isDrift,0.);
2650 INTERP_KERNEL::matrixProduct(matrixInv->getConstPointer(),nbOfPts+isDrift,nbOfPts+isDrift,arr2->getConstPointer(),nbOfPts+isDrift,1,KnewiK->getPointer());
2651 return KnewiK.retn();
2655 * 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.
2657 * \param [in] spaceDimension space dimension of the input mesh on which the Kriging has to be performed
2658 * \param [in] nbOfElems is the result of the product of nb of rows and the nb of columns of matrix \a matrixPtr
2659 * \param [in,out] matrixPtr is the dense matrix whose on each values the operation will be applied
2661 void MEDCouplingFieldDiscretizationKriging::operateOnDenseMatrix(int spaceDimension, int nbOfElems, double *matrixPtr) const
2663 switch(spaceDimension)
2667 for(int i=0;i<nbOfElems;i++)
2669 double val=matrixPtr[i];
2670 matrixPtr[i]=val*val*val;
2675 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationKriging::operateOnDenseMatrix : only dimension 1 implemented !");
2680 * Starting from a square matrix \a matr, this method returns a newly allocated dense square matrix whose \a matr is included in returned matrix
2681 * in the top left corner, and in the remaining returned matrix the parameters to take into account about the kriging drift.
2682 * For the moment only linear srift is implemented.
2684 * \param [in] arr the position of points were input mesh geometry is considered for Kriging
2685 * \param [in] matr input matrix whose drift part will be added
2686 * \param [out] delta the difference between the size of the output matrix and the input matrix \a matr.
2687 * \return a newly allocated matrix bigger than input matrix \a matr.
2689 DataArrayDouble *MEDCouplingFieldDiscretizationKriging::performDrift(const DataArrayDouble *matr, const DataArrayDouble *arr, int& delta) const
2691 int spaceDimension=arr->getNumberOfComponents();
2692 delta=spaceDimension+1;
2693 int szOfMatrix=arr->getNumberOfTuples();
2694 if(szOfMatrix*szOfMatrix!=matr->getNumberOfTuples())
2695 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationKriging::performDrift : invalid size");
2696 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
2697 ret->alloc((szOfMatrix+delta)*(szOfMatrix+delta),1);
2698 const double *srcWork=matr->getConstPointer();
2699 const double *srcWork2=arr->getConstPointer();
2700 double *destWork=ret->getPointer();
2701 for(int i=0;i<szOfMatrix;i++)
2703 destWork=std::copy(srcWork,srcWork+szOfMatrix,destWork);
2704 srcWork+=szOfMatrix;
2706 destWork=std::copy(srcWork2,srcWork2+spaceDimension,destWork);
2707 srcWork2+=spaceDimension;
2709 std::fill(destWork,destWork+szOfMatrix,1.); destWork+=szOfMatrix;
2710 std::fill(destWork,destWork+spaceDimension+1,0.); destWork+=spaceDimension+1;
2711 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arrNoI=arr->toNoInterlace();
2712 srcWork2=arrNoI->getConstPointer();
2713 for(int i=0;i<spaceDimension;i++)
2715 destWork=std::copy(srcWork2,srcWork2+szOfMatrix,destWork);
2716 srcWork2+=szOfMatrix;
2717 std::fill(destWork,destWork+spaceDimension+1,0.);
2718 destWork+=spaceDimension;