1 // Copyright (C) 2007-2014 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, or (at your option) any later version.
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.8888888888888888,0.5555555555555556};
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_QUAD8[8]={1.,1.,1.,1.,1.,1.,1.,1.};
76 const double MEDCouplingFieldDiscretizationGaussNE::FGP_QUAD9[9]={0.30864197530864196,0.30864197530864196,0.30864197530864196,0.30864197530864196,0.49382716049382713,0.49382716049382713,0.49382716049382713,0.49382716049382713,0.7901234567901234};
77 const double MEDCouplingFieldDiscretizationGaussNE::FGP_TETRA4[4]={0.041666666666666664,0.041666666666666664,0.041666666666666664,0.041666666666666664};
78 const double MEDCouplingFieldDiscretizationGaussNE::FGP_PENTA6[6]={0.16666666666666666,0.16666666666666666,0.16666666666666666,0.16666666666666666,0.16666666666666666,0.16666666666666666};
79 const double MEDCouplingFieldDiscretizationGaussNE::FGP_HEXA8[8]={1.,1.,1.,1.,1.,1.,1.,1.};
80 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};
81 const double MEDCouplingFieldDiscretizationGaussNE::FGP_PYRA5[5]={0.13333333333333333,0.13333333333333333,0.13333333333333333,0.13333333333333333,0.13333333333333333};
82 const double MEDCouplingFieldDiscretizationGaussNE::REF_SEG2[2]={-1.,1.};
83 const double MEDCouplingFieldDiscretizationGaussNE::REF_SEG3[3]={-1.,1.,0.};
84 const double MEDCouplingFieldDiscretizationGaussNE::REF_SEG4[4]={-1.,1.,-0.3333333333333333,0.3333333333333333};
85 const double MEDCouplingFieldDiscretizationGaussNE::REF_TRI3[6]={0.,0.,1.,0.,0.,1.};
86 const double MEDCouplingFieldDiscretizationGaussNE::REF_TRI6[12]={0.,0.,1.,0.,0.,1.,0.5,0.,0.5,0.5,0.,0.5};
87 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};
88 const double MEDCouplingFieldDiscretizationGaussNE::REF_QUAD4[8]={-1.,-1.,1.,-1.,1.,1.,-1.,1.};
89 const double MEDCouplingFieldDiscretizationGaussNE::REF_QUAD8[16]={-1.,-1.,1.,-1.,1.,1.,-1.,1.,0.,-1.,1.,0.,0.,1.,-1.,0.};
90 const double MEDCouplingFieldDiscretizationGaussNE::REF_QUAD9[18]={-1.,-1.,1.,-1.,1.,1.,-1.,1.,0.,-1.,1.,0.,0.,1.,-1.,0.,0.,0.};
91 const double MEDCouplingFieldDiscretizationGaussNE::REF_TETRA4[12]={0.,1.,0.,0.,0.,1.,0.,0.,0.,1.,0.,0.};
92 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.};
93 const double MEDCouplingFieldDiscretizationGaussNE::REF_PENTA6[18]={-1.,1.,0.,-1.,0.,1.,-1.,0.,0.,1.,1.,0.,1.,0.,1.,1.,0.,0.};
94 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.};
95 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.};
96 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.};
97 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.};
98 const double MEDCouplingFieldDiscretizationGaussNE::REF_PYRA5[15]={1.,0.,0.,0.,1.,0.,-1.,0.,0.,0.,-1.,0.,0.,0.,1.};
99 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 const double MEDCouplingFieldDiscretizationGaussNE::LOC_SEG2[2]={0.577350269189626,-0.577350269189626};
101 const double MEDCouplingFieldDiscretizationGaussNE::LOC_SEG3[3]={-0.774596669241,0.,0.774596669241};
102 const double MEDCouplingFieldDiscretizationGaussNE::LOC_SEG4[4]={0.339981043584856,-0.339981043584856,0.861136311594053,-0.861136311594053};
103 const double MEDCouplingFieldDiscretizationGaussNE::LOC_TRI3[6]={0.16666666666666667,0.16666666666666667,0.6666666666666667,0.16666666666666667,0.16666666666666667,0.6666666666666667};
104 const double MEDCouplingFieldDiscretizationGaussNE::LOC_TRI6[12]={0.091576213509771,0.091576213509771,0.816847572980458,0.091576213509771,0.091576213509771,0.816847572980458,0.445948490915965,0.10810301816807,0.445948490915965,0.445948490915965,0.10810301816807,0.445948490915965};
105 const double MEDCouplingFieldDiscretizationGaussNE::LOC_TRI7[14]={0.3333333333333333,0.3333333333333333,0.470142064105115,0.470142064105115,0.05971587178977,0.470142064105115,0.470142064105115,0.05971587178977,0.101286507323456,0.101286507323456,0.797426985353088,0.101286507323456,0.101286507323456,0.797426985353088};
106 const double MEDCouplingFieldDiscretizationGaussNE::LOC_QUAD4[8]={-0.774596669241483,-0.774596669241483,0.774596669241483,-0.774596669241483,0.774596669241483,0.774596669241483,-0.774596669241483,0.774596669241483};
107 const double MEDCouplingFieldDiscretizationGaussNE::LOC_QUAD8[16]={-0.774596669241483,-0.774596669241483,0.774596669241483,-0.774596669241483,0.774596669241483,0.774596669241483,-0.774596669241483,0.774596669241483,0.,-0.774596669241483,0.774596669241483,0.,0.,0.774596669241483,-0.774596669241483,0.};
108 const double MEDCouplingFieldDiscretizationGaussNE::LOC_QUAD9[18]={-0.774596669241483,-0.774596669241483,0.774596669241483,-0.774596669241483,0.774596669241483,0.774596669241483,-0.774596669241483,0.774596669241483,0.,-0.774596669241483,0.774596669241483,0.,0.,0.774596669241483,-0.774596669241483,0.,0.,0.};
109 const double MEDCouplingFieldDiscretizationGaussNE::LOC_TETRA4[12]={0.1381966011250105,0.1381966011250105,0.1381966011250105,0.1381966011250105,0.1381966011250105,0.5854101966249685,0.1381966011250105,0.5854101966249685,0.1381966011250105,0.5854101966249685,0.1381966011250105,0.1381966011250105};
110 const double MEDCouplingFieldDiscretizationGaussNE::LOC_PENTA6[18]={-0.5773502691896258,0.5,0.5,-0.5773502691896258,0.,0.5,-0.5773502691896258,0.5,0.,0.5773502691896258,0.5,0.5,0.5773502691896258,0.,0.5,0.5773502691896258,0.5,0.};
111 const double MEDCouplingFieldDiscretizationGaussNE::LOC_HEXA8[24]={-0.5773502691896258,-0.5773502691896258,-0.5773502691896258,-0.5773502691896258,-0.5773502691896258,0.5773502691896258,-0.5773502691896258,0.5773502691896258,-0.5773502691896258,-0.5773502691896258,0.5773502691896258,0.5773502691896258,0.5773502691896258,-0.5773502691896258,-0.5773502691896258,0.5773502691896258,-0.5773502691896258,0.5773502691896258,0.5773502691896258,0.5773502691896258,-0.5773502691896258,0.5773502691896258,0.5773502691896258,0.5773502691896258};
112 const double MEDCouplingFieldDiscretizationGaussNE::LOC_HEXA27[81]={-0.7745966692414834,-0.7745966692414834,-0.7745966692414834,-0.7745966692414834,-0.7745966692414834,0.,-0.7745966692414834,-0.7745966692414834,0.7745966692414834,-0.7745966692414834,0.,-0.7745966692414834,-0.7745966692414834,0.,0.,-0.7745966692414834,0.,0.7745966692414834,-0.7745966692414834,0.7745966692414834,-0.7745966692414834,-0.7745966692414834,0.7745966692414834,0.,-0.7745966692414834,0.7745966692414834,0.7745966692414834,0.,-0.7745966692414834,-0.7745966692414834,0,-0.7745966692414834,0.,0.,-0.7745966692414834,0.7745966692414834,0.,0.,-0.7745966692414834,0.,0.,0.,0.,0.,0.7745966692414834,0.,0.7745966692414834,-0.7745966692414834,0.,0.7745966692414834,0.,0.,0.7745966692414834,0.7745966692414834,0.7745966692414834,-0.7745966692414834,-0.7745966692414834,0.7745966692414834,-0.7745966692414834,0.,0.7745966692414834,-0.7745966692414834,0.7745966692414834,0.7745966692414834,0,-0.7745966692414834,0.7745966692414834,0.,0.,0.7745966692414834,0.,0.7745966692414834,0.7745966692414834,0.7745966692414834,-0.7745966692414834,0.7745966692414834,0.7745966692414834,0.,0.7745966692414834,0.7745966692414834,0.7745966692414834};
113 const double MEDCouplingFieldDiscretizationGaussNE::LOC_PYRA5[15]={0.5,0.,0.1531754163448146,0.,0.5,0.1531754163448146,-0.5,0.,0.1531754163448146,0.,-0.5,0.1531754163448146,0.,0.,0.6372983346207416};
115 MEDCouplingFieldDiscretization::MEDCouplingFieldDiscretization():_precision(DFLT_PRECISION)
119 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretization::New(TypeOfField type)
123 case MEDCouplingFieldDiscretizationP0::TYPE:
124 return new MEDCouplingFieldDiscretizationP0;
125 case MEDCouplingFieldDiscretizationP1::TYPE:
126 return new MEDCouplingFieldDiscretizationP1;
127 case MEDCouplingFieldDiscretizationGauss::TYPE:
128 return new MEDCouplingFieldDiscretizationGauss;
129 case MEDCouplingFieldDiscretizationGaussNE::TYPE:
130 return new MEDCouplingFieldDiscretizationGaussNE;
131 case MEDCouplingFieldDiscretizationKriging::TYPE:
132 return new MEDCouplingFieldDiscretizationKriging;
134 throw INTERP_KERNEL::Exception("Choosen discretization is not implemented yet.");
138 TypeOfField MEDCouplingFieldDiscretization::GetTypeOfFieldFromStringRepr(const std::string& repr)
140 if(repr==MEDCouplingFieldDiscretizationP0::REPR)
141 return MEDCouplingFieldDiscretizationP0::TYPE;
142 if(repr==MEDCouplingFieldDiscretizationP1::REPR)
143 return MEDCouplingFieldDiscretizationP1::TYPE;
144 if(repr==MEDCouplingFieldDiscretizationGauss::REPR)
145 return MEDCouplingFieldDiscretizationGauss::TYPE;
146 if(repr==MEDCouplingFieldDiscretizationGaussNE::REPR)
147 return MEDCouplingFieldDiscretizationGaussNE::TYPE;
148 if(repr==MEDCouplingFieldDiscretizationKriging::REPR)
149 return MEDCouplingFieldDiscretizationKriging::TYPE;
150 throw INTERP_KERNEL::Exception("Representation does not match with any field discretization !");
153 bool MEDCouplingFieldDiscretization::isEqual(const MEDCouplingFieldDiscretization *other, double eps) const
156 return isEqualIfNotWhy(other,eps,reason);
159 bool MEDCouplingFieldDiscretization::isEqualWithoutConsideringStr(const MEDCouplingFieldDiscretization *other, double eps) const
161 return isEqual(other,eps);
165 * This method is an alias of MEDCouplingFieldDiscretization::clone. It is only here for coherency with all the remaining of MEDCoupling.
166 * \sa MEDCouplingFieldDiscretization::clone.
168 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretization::deepCpy() const
174 * For all field discretization excepted GaussPts the [ \a startCellIds, \a endCellIds ) has no impact on the cloned instance.
176 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretization::clonePart(const int *startCellIds, const int *endCellIds) const
182 * For all field discretization excepted GaussPts the slice( \a beginCellId, \a endCellIds, \a stepCellId ) has no impact on the cloned instance.
184 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretization::clonePartRange(int beginCellIds, int endCellIds, int stepCellIds) const
190 * Excepted for MEDCouplingFieldDiscretizationPerCell no underlying TimeLabel object : nothing to do in generally.
192 void MEDCouplingFieldDiscretization::updateTime() const
196 std::size_t MEDCouplingFieldDiscretization::getHeapMemorySizeWithoutChildren() const
201 std::vector<const BigMemoryObject *> MEDCouplingFieldDiscretization::getDirectChildren() const
203 return std::vector<const BigMemoryObject *>();
207 * Computes normL1 of DataArrayDouble instance arr.
208 * @param res output parameter expected to be of size arr->getNumberOfComponents();
209 * @throw when the field discretization fails on getMeasure fields (gauss points for example)
211 void MEDCouplingFieldDiscretization::normL1(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, double *res) const
213 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=getMeasureField(mesh,true);
214 int nbOfCompo=arr->getNumberOfComponents();
215 int nbOfElems=getNumberOfTuples(mesh);
216 std::fill(res,res+nbOfCompo,0.);
217 const double *arrPtr=arr->getConstPointer();
218 const double *volPtr=vol->getArray()->getConstPointer();
220 for(int i=0;i<nbOfElems;i++)
222 double v=fabs(volPtr[i]);
223 for(int j=0;j<nbOfCompo;j++)
224 res[j]+=fabs(arrPtr[i*nbOfCompo+j])*v;
227 std::transform(res,res+nbOfCompo,res,std::bind2nd(std::multiplies<double>(),1./deno));
231 * Computes normL2 of DataArrayDouble instance arr.
232 * @param res output parameter expected to be of size arr->getNumberOfComponents();
233 * @throw when the field discretization fails on getMeasure fields (gauss points for example)
235 void MEDCouplingFieldDiscretization::normL2(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, double *res) const
237 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=getMeasureField(mesh,true);
238 int nbOfCompo=arr->getNumberOfComponents();
239 int nbOfElems=getNumberOfTuples(mesh);
240 std::fill(res,res+nbOfCompo,0.);
241 const double *arrPtr=arr->getConstPointer();
242 const double *volPtr=vol->getArray()->getConstPointer();
244 for(int i=0;i<nbOfElems;i++)
246 double v=fabs(volPtr[i]);
247 for(int j=0;j<nbOfCompo;j++)
248 res[j]+=arrPtr[i*nbOfCompo+j]*arrPtr[i*nbOfCompo+j]*v;
251 std::transform(res,res+nbOfCompo,res,std::bind2nd(std::multiplies<double>(),1./deno));
252 std::transform(res,res+nbOfCompo,res,std::ptr_fun<double,double>(std::sqrt));
256 * Computes integral of DataArrayDouble instance arr.
257 * @param res output parameter expected to be of size arr->getNumberOfComponents();
258 * @throw when the field discretization fails on getMeasure fields (gauss points for example)
260 void MEDCouplingFieldDiscretization::integral(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, bool isWAbs, double *res) const
263 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretization::integral : mesh is NULL !");
265 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretization::integral : input array is NULL !");
266 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=getMeasureField(mesh,isWAbs);
267 int nbOfCompo=arr->getNumberOfComponents();
268 int nbOfElems=getNumberOfTuples(mesh);
269 if(nbOfElems!=arr->getNumberOfTuples())
271 std::ostringstream oss; oss << "MEDCouplingFieldDiscretization::integral : field is not correct ! number of tuples in array is " << arr->getNumberOfTuples();
272 oss << " whereas number of tuples expected is " << nbOfElems << " !";
273 throw INTERP_KERNEL::Exception(oss.str().c_str());
275 std::fill(res,res+nbOfCompo,0.);
276 const double *arrPtr=arr->getConstPointer();
277 const double *volPtr=vol->getArray()->getConstPointer();
278 INTERP_KERNEL::AutoPtr<double> tmp=new double[nbOfCompo];
279 for (int i=0;i<nbOfElems;i++)
281 std::transform(arrPtr+i*nbOfCompo,arrPtr+(i+1)*nbOfCompo,(double *)tmp,std::bind2nd(std::multiplies<double>(),volPtr[i]));
282 std::transform((double *)tmp,(double *)tmp+nbOfCompo,res,res,std::plus<double>());
287 * This method is strictly equivalent to MEDCouplingFieldDiscretization::buildSubMeshData except that it is optimized for input defined as a range of cell ids.
289 * \param [out] beginOut Valid only if \a di is NULL
290 * \param [out] endOut Valid only if \a di is NULL
291 * \param [out] stepOut Valid only if \a di is NULL
292 * \param [out] di is an array returned that specifies entity ids (nodes, cells, Gauss points... ) in array if no output range is foundable.
294 * \sa MEDCouplingFieldDiscretization::buildSubMeshData
296 MEDCouplingMesh *MEDCouplingFieldDiscretization::buildSubMeshDataRange(const MEDCouplingMesh *mesh, int beginCellIds, int endCellIds, int stepCellIds, int& beginOut, int& endOut, int& stepOut, DataArrayInt *&di) const
298 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> da=DataArrayInt::Range(beginCellIds,endCellIds,stepCellIds);
299 return buildSubMeshData(mesh,da->begin(),da->end(),di);
302 void MEDCouplingFieldDiscretization::getSerializationIntArray(DataArrayInt *& arr) const
310 void MEDCouplingFieldDiscretization::getTinySerializationIntInformation(std::vector<int>& tinyInfo) const
317 void MEDCouplingFieldDiscretization::getTinySerializationDbleInformation(std::vector<double>& tinyInfo) const
321 void MEDCouplingFieldDiscretization::resizeForUnserialization(const std::vector<int>& tinyInfo, DataArrayInt *& arr)
329 void MEDCouplingFieldDiscretization::finishUnserialization(const std::vector<double>& tinyInfo)
334 * This method is typically the first step of renumbering. The implementation is empty it is not a bug only gauss is impacted
335 * virtualy by this method.
337 void MEDCouplingFieldDiscretization::renumberCells(const int *old2NewBg, bool check)
341 double MEDCouplingFieldDiscretization::getIJK(const MEDCouplingMesh *mesh, const DataArrayDouble *da, int cellId, int nodeIdInCell, int compoId) const
343 throw INTERP_KERNEL::Exception("getIJK Invalid ! only for GaussPoint and GaussNE discretizations !");
346 void MEDCouplingFieldDiscretization::setGaussLocalizationOnType(const MEDCouplingMesh *m, INTERP_KERNEL::NormalizedCellType type, const std::vector<double>& refCoo,
347 const std::vector<double>& gsCoo, const std::vector<double>& wg) throw(INTERP_KERNEL::Exception)
349 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
352 void MEDCouplingFieldDiscretization::setGaussLocalizationOnCells(const MEDCouplingMesh *m, const int *begin, const int *end, const std::vector<double>& refCoo,
353 const std::vector<double>& gsCoo, const std::vector<double>& wg) throw(INTERP_KERNEL::Exception)
355 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
358 void MEDCouplingFieldDiscretization::clearGaussLocalizations()
360 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
363 MEDCouplingGaussLocalization& MEDCouplingFieldDiscretization::getGaussLocalization(int locId)
365 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
368 const MEDCouplingGaussLocalization& MEDCouplingFieldDiscretization::getGaussLocalization(int locId) const
370 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
373 int MEDCouplingFieldDiscretization::getNbOfGaussLocalization() const
375 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
378 int MEDCouplingFieldDiscretization::getGaussLocalizationIdOfOneCell(int cellId) const
380 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
383 int MEDCouplingFieldDiscretization::getGaussLocalizationIdOfOneType(INTERP_KERNEL::NormalizedCellType type) const
385 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
388 std::set<int> MEDCouplingFieldDiscretization::getGaussLocalizationIdsOfOneType(INTERP_KERNEL::NormalizedCellType type) const
390 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
393 void MEDCouplingFieldDiscretization::getCellIdsHavingGaussLocalization(int locId, std::vector<int>& cellIds) const
395 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
398 void MEDCouplingFieldDiscretization::RenumberEntitiesFromO2NArr(double eps, const int *old2NewPtr, int newNbOfEntity, DataArrayDouble *arr, const std::string& msg)
401 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretization::RenumberEntitiesFromO2NArr : input array is NULL !");
402 int oldNbOfElems=arr->getNumberOfTuples();
403 int nbOfComp=arr->getNumberOfComponents();
404 int newNbOfTuples=newNbOfEntity;
405 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arrCpy=arr->deepCpy();
406 const double *ptSrc=arrCpy->getConstPointer();
407 arr->reAlloc(newNbOfTuples);
408 double *ptToFill=arr->getPointer();
409 std::fill(ptToFill,ptToFill+nbOfComp*newNbOfTuples,std::numeric_limits<double>::max());
410 INTERP_KERNEL::AutoPtr<double> tmp=new double[nbOfComp];
411 for(int i=0;i<oldNbOfElems;i++)
413 int newNb=old2NewPtr[i];
414 if(newNb>=0)//if newNb<0 the node is considered as out.
416 if(std::find_if(ptToFill+newNb*nbOfComp,ptToFill+(newNb+1)*nbOfComp,std::bind2nd(std::not_equal_to<double>(),std::numeric_limits<double>::max()))
417 ==ptToFill+(newNb+1)*nbOfComp)
418 std::copy(ptSrc+i*nbOfComp,ptSrc+(i+1)*nbOfComp,ptToFill+newNb*nbOfComp);
421 std::transform(ptSrc+i*nbOfComp,ptSrc+(i+1)*nbOfComp,ptToFill+newNb*nbOfComp,(double *)tmp,std::minus<double>());
422 std::transform((double *)tmp,((double *)tmp)+nbOfComp,(double *)tmp,std::ptr_fun<double,double>(fabs));
423 //if(!std::equal(ptSrc+i*nbOfComp,ptSrc+(i+1)*nbOfComp,ptToFill+newNb*nbOfComp))
424 if(*std::max_element((double *)tmp,((double *)tmp)+nbOfComp)>eps)
426 std::ostringstream oss;
427 oss << msg << " " << i << " and " << std::find(old2NewPtr,old2NewPtr+i,newNb)-old2NewPtr
428 << " have been merged and " << msg << " field on them are different !";
429 throw INTERP_KERNEL::Exception(oss.str().c_str());
436 void MEDCouplingFieldDiscretization::RenumberEntitiesFromN2OArr(const int *new2OldPtr, int new2OldSz, DataArrayDouble *arr, const std::string& msg)
438 int nbOfComp=arr->getNumberOfComponents();
439 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arrCpy=arr->deepCpy();
440 const double *ptSrc=arrCpy->getConstPointer();
441 arr->reAlloc(new2OldSz);
442 double *ptToFill=arr->getPointer();
443 for(int i=0;i<new2OldSz;i++)
445 int oldNb=new2OldPtr[i];
446 std::copy(ptSrc+oldNb*nbOfComp,ptSrc+(oldNb+1)*nbOfComp,ptToFill+i*nbOfComp);
450 MEDCouplingFieldDiscretization::~MEDCouplingFieldDiscretization()
454 TypeOfField MEDCouplingFieldDiscretizationP0::getEnum() const
460 * This method is simply called by MEDCouplingFieldDiscretization::deepCpy. It performs the deep copy of \a this.
462 * \sa MEDCouplingFieldDiscretization::deepCpy.
464 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationP0::clone() const
466 return new MEDCouplingFieldDiscretizationP0;
469 std::string MEDCouplingFieldDiscretizationP0::getStringRepr() const
471 return std::string(REPR);
474 const char *MEDCouplingFieldDiscretizationP0::getRepr() const
479 bool MEDCouplingFieldDiscretizationP0::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
483 reason="other spatial discretization is NULL, and this spatial discretization (P0) is defined.";
486 const MEDCouplingFieldDiscretizationP0 *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationP0 *>(other);
489 reason="Spatial discrtization of this is ON_CELLS, which is not the case of other.";
493 int MEDCouplingFieldDiscretizationP0::getNumberOfTuples(const MEDCouplingMesh *mesh) const
496 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getNumberOfTuples : NULL input mesh !");
497 return mesh->getNumberOfCells();
501 * This method returns the number of tuples regarding exclusively the input code \b without \b using \b a \b mesh \b in \b input.
502 * The input code coherency is also checked regarding spatial discretization of \a this.
503 * If an incoherency is detected, an exception will be thrown. If the input code is coherent, the number of tuples expected is returned.
504 * The number of tuples expected is equal to those to have a valid field lying on \a this and having a mesh fitting perfectly the input code (geometric type distribution).
506 int MEDCouplingFieldDiscretizationP0::getNumberOfTuplesExpectedRegardingCode(const std::vector<int>& code, const std::vector<const DataArrayInt *>& idsPerType) const
509 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getNumberOfTuplesExpectedRegardingCode : invalid input code !");
510 int nbOfSplit=(int)idsPerType.size();
511 int nbOfTypes=(int)code.size()/3;
513 for(int i=0;i<nbOfTypes;i++)
515 int nbOfEltInChunk=code[3*i+1];
517 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getNumberOfTuplesExpectedRegardingCode : invalid input code ! presence of negative value in a type !");
521 if(pos<0 || pos>=nbOfSplit)
523 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationP0::getNumberOfTuplesExpectedRegardingCode : input code points to pos " << pos << " in typeid " << i << " ! Should be in [0," << nbOfSplit << ") !";
524 throw INTERP_KERNEL::Exception(oss.str().c_str());
526 const DataArrayInt *ids(idsPerType[pos]);
527 if(!ids || !ids->isAllocated() || ids->getNumberOfComponents()!=1 || ids->getNumberOfTuples()!=nbOfEltInChunk || ids->getMinValueInArray()<0)
529 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationP0::getNumberOfTuplesExpectedRegardingCode : input pfl chunck at pos " << pos << " should have " << i << " tuples and one component and with ids all >=0 !";
530 throw INTERP_KERNEL::Exception(oss.str().c_str());
538 int MEDCouplingFieldDiscretizationP0::getNumberOfMeshPlaces(const MEDCouplingMesh *mesh) const
541 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getNumberOfMeshPlaces : NULL input mesh !");
542 return mesh->getNumberOfCells();
545 DataArrayInt *MEDCouplingFieldDiscretizationP0::getOffsetArr(const MEDCouplingMesh *mesh) const
548 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getOffsetArr : NULL input mesh !");
549 int nbOfTuples=mesh->getNumberOfCells();
550 DataArrayInt *ret=DataArrayInt::New();
551 ret->alloc(nbOfTuples+1,1);
556 void MEDCouplingFieldDiscretizationP0::renumberArraysForCell(const MEDCouplingMesh *mesh, const std::vector<DataArray *>& arrays,
557 const int *old2NewBg, bool check) throw(INTERP_KERNEL::Exception)
560 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::renumberArraysForCell : NULL input mesh !");
561 const int *array=old2NewBg;
563 array=DataArrayInt::CheckAndPreparePermutation(old2NewBg,old2NewBg+mesh->getNumberOfCells());
564 for(std::vector<DataArray *>::const_iterator it=arrays.begin();it!=arrays.end();it++)
567 (*it)->renumberInPlace(array);
570 free(const_cast<int *>(array));
573 DataArrayDouble *MEDCouplingFieldDiscretizationP0::getLocalizationOfDiscValues(const MEDCouplingMesh *mesh) const
576 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getLocalizationOfDiscValues : NULL input mesh !");
577 return mesh->getBarycenterAndOwner();
580 void MEDCouplingFieldDiscretizationP0::computeMeshRestrictionFromTupleIds(const MEDCouplingMesh *mesh, const int *tupleIdsBg, const int *tupleIdsEnd,
581 DataArrayInt *&cellRestriction, DataArrayInt *&trueTupleRestriction) const throw(INTERP_KERNEL::Exception)
584 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::computeMeshRestrictionFromTupleIds : NULL input mesh !");
585 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> tmp=DataArrayInt::New();
586 tmp->alloc((int)std::distance(tupleIdsBg,tupleIdsEnd),1);
587 std::copy(tupleIdsBg,tupleIdsEnd,tmp->getPointer());
588 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> tmp2(tmp->deepCpy());
589 cellRestriction=tmp.retn();
590 trueTupleRestriction=tmp2.retn();
593 void MEDCouplingFieldDiscretizationP0::reprQuickOverview(std::ostream& stream) const
595 stream << "P0 spatial discretization.";
598 void MEDCouplingFieldDiscretizationP0::checkCompatibilityWithNature(NatureOfField nat) const
602 void MEDCouplingFieldDiscretizationP0::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArray *da) const
605 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::checkCoherencyBetween : NULL input mesh or DataArray !");
606 if(mesh->getNumberOfCells()!=da->getNumberOfTuples())
608 std::ostringstream message;
609 message << "Field on cells invalid because there are " << mesh->getNumberOfCells();
610 message << " cells in mesh and " << da->getNumberOfTuples() << " tuples in field !";
611 throw INTERP_KERNEL::Exception(message.str().c_str());
615 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationP0::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
618 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getMeasureField : mesh instance specified is NULL !");
619 return mesh->getMeasureField(isAbs);
622 void MEDCouplingFieldDiscretizationP0::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
625 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getValueOn : NULL input mesh !");
626 int id=mesh->getCellContainingPoint(loc,_precision);
628 throw INTERP_KERNEL::Exception("Specified point is detected outside of mesh : unable to apply P0::getValueOn !");
629 arr->getTuple(id,res);
632 void MEDCouplingFieldDiscretizationP0::getValueOnPos(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, int i, int j, int k, double *res) const
634 const MEDCouplingCMesh *meshC=dynamic_cast<const MEDCouplingCMesh *>(mesh);
636 throw INTERP_KERNEL::Exception("P0::getValueOnPos is only accessible for structured meshes !");
637 int id=meshC->getCellIdFromPos(i,j,k);
638 arr->getTuple(id,res);
641 DataArrayDouble *MEDCouplingFieldDiscretizationP0::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const
644 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getValueOnMulti : NULL input mesh !");
645 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> eltsArr,eltsIndexArr;
646 mesh->getCellsContainingPoints(loc,nbOfPoints,_precision,eltsArr,eltsIndexArr);
647 const int *elts(eltsArr->begin()),*eltsIndex(eltsIndexArr->begin());
648 int spaceDim=mesh->getSpaceDimension();
649 int nbOfComponents=arr->getNumberOfComponents();
650 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
651 ret->alloc(nbOfPoints,nbOfComponents);
652 double *ptToFill=ret->getPointer();
653 for(int i=0;i<nbOfPoints;i++,ptToFill+=nbOfComponents)
654 if(eltsIndex[i+1]-eltsIndex[i]>=1)
655 arr->getTuple(elts[eltsIndex[i]],ptToFill);
658 std::ostringstream oss; oss << "Point #" << i << " with coordinates : (";
659 std::copy(loc+i*spaceDim,loc+(i+1)*spaceDim,std::ostream_iterator<double>(oss,", "));
660 oss << ") detected outside mesh : unable to apply P0::getValueOnMulti ! ";
661 throw INTERP_KERNEL::Exception(oss.str().c_str());
667 * Nothing to do. It's not a bug.
669 void MEDCouplingFieldDiscretizationP0::renumberValuesOnNodes(double , const int *, int newNbOfNodes, DataArrayDouble *) const
673 void MEDCouplingFieldDiscretizationP0::renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const
675 RenumberEntitiesFromO2NArr(epsOnVals,old2New,newSz,arr,"Cell");
678 void MEDCouplingFieldDiscretizationP0::renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const
680 RenumberEntitiesFromN2OArr(new2old,newSz,arr,"Cell");
684 * This method returns a tuple ids selection from cell ids selection [start;end).
685 * This method is called by MEDCouplingFieldDiscretizationP0::buildSubMeshData to return parameter \b di.
686 * Here for P0 it's very simple !
688 * \return a newly allocated array containing ids to select into the DataArrayDouble of the field.
691 DataArrayInt *MEDCouplingFieldDiscretizationP0::computeTupleIdsToSelectFromCellIds(const MEDCouplingMesh *mesh, const int *startCellIds, const int *endCellIds) const
693 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret=DataArrayInt::New();
694 ret->alloc((int)std::distance(startCellIds,endCellIds),1);
695 std::copy(startCellIds,endCellIds,ret->getPointer());
700 * This method returns a submesh of 'mesh' instance constituting cell ids contained in array defined as an interval [start;end).
701 * @param di is an array returned that specifies entity ids (here cells ids) in mesh 'mesh' of entity in returned submesh.
702 * Example : The first cell id of returned mesh has the (*di)[0] id in 'mesh'
704 * \sa MEDCouplingFieldDiscretizationP0::buildSubMeshDataRange
706 MEDCouplingMesh *MEDCouplingFieldDiscretizationP0::buildSubMeshData(const MEDCouplingMesh *mesh, const int *start, const int *end, DataArrayInt *&di) const
709 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::buildSubMeshData : NULL input mesh !");
710 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPart(start,end);
711 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> diSafe=DataArrayInt::New();
712 diSafe->alloc((int)std::distance(start,end),1);
713 std::copy(start,end,diSafe->getPointer());
719 * This method is strictly equivalent to MEDCouplingFieldDiscretizationP0::buildSubMeshData except that it is optimized for input defined as a range of cell ids.
721 * \param [out] beginOut Valid only if \a di is NULL
722 * \param [out] endOut Valid only if \a di is NULL
723 * \param [out] stepOut Valid only if \a di is NULL
724 * \param [out] di is an array returned that specifies entity ids (nodes, cells, Gauss points... ) in array if no output range is foundable.
726 * \sa MEDCouplingFieldDiscretizationP0::buildSubMeshData
728 MEDCouplingMesh *MEDCouplingFieldDiscretizationP0::buildSubMeshDataRange(const MEDCouplingMesh *mesh, int beginCellIds, int endCellIds, int stepCellIds, int& beginOut, int& endOut, int& stepOut, DataArrayInt *&di) const
731 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::buildSubMeshDataRange : NULL input mesh !");
732 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPartRange(beginCellIds,endCellIds,stepCellIds);
733 di=0; beginOut=beginCellIds; endOut=endCellIds; stepOut=stepCellIds;
737 int MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuples(const MEDCouplingMesh *mesh) const
740 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::getNumberOfTuples : NULL input mesh !");
741 return mesh->getNumberOfNodes();
745 * This method returns the number of tuples regarding exclusively the input code \b without \b using \b a \b mesh \b in \b input.
746 * The input code coherency is also checked regarding spatial discretization of \a this.
747 * If an incoherency is detected, an exception will be thrown. If the input code is coherent, the number of tuples expected is returned.
748 * The number of tuples expected is equal to those to have a valid field lying on \a this and having a mesh fitting perfectly the input code (geometric type distribution).
750 int MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuplesExpectedRegardingCode(const std::vector<int>& code, const std::vector<const DataArrayInt *>& idsPerType) const
753 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuplesExpectedRegardingCode : invalid input code !");
754 int nbOfSplit=(int)idsPerType.size();
755 int nbOfTypes=(int)code.size()/3;
757 for(int i=0;i<nbOfTypes;i++)
759 int nbOfEltInChunk=code[3*i+1];
761 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuplesExpectedRegardingCode : invalid input code ! presence of negative value in a type !");
765 if(pos<0 || pos>=nbOfSplit)
767 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuplesExpectedRegardingCode : input code points to pos " << pos << " in typeid " << i << " ! Should be in [0," << nbOfSplit << ") !";
768 throw INTERP_KERNEL::Exception(oss.str().c_str());
770 const DataArrayInt *ids(idsPerType[pos]);
771 if(!ids || !ids->isAllocated() || ids->getNumberOfComponents()!=1 || ids->getNumberOfTuples()!=nbOfEltInChunk || ids->getMinValueInArray()<0)
773 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuplesExpectedRegardingCode : input pfl chunck at pos " << pos << " should have " << i << " tuples and one component and with ids all >=0 !";
774 throw INTERP_KERNEL::Exception(oss.str().c_str());
782 int MEDCouplingFieldDiscretizationOnNodes::getNumberOfMeshPlaces(const MEDCouplingMesh *mesh) const
785 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::getNumberOfMeshPlaces : NULL input mesh !");
786 return mesh->getNumberOfNodes();
790 * Nothing to do here.
792 void MEDCouplingFieldDiscretizationOnNodes::renumberArraysForCell(const MEDCouplingMesh *, const std::vector<DataArray *>& arrays,
793 const int *old2NewBg, bool check) throw(INTERP_KERNEL::Exception)
797 DataArrayInt *MEDCouplingFieldDiscretizationOnNodes::getOffsetArr(const MEDCouplingMesh *mesh) const
800 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::getOffsetArr : NULL input mesh !");
801 int nbOfTuples=mesh->getNumberOfNodes();
802 DataArrayInt *ret=DataArrayInt::New();
803 ret->alloc(nbOfTuples+1,1);
808 DataArrayDouble *MEDCouplingFieldDiscretizationOnNodes::getLocalizationOfDiscValues(const MEDCouplingMesh *mesh) const
811 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::getLocalizationOfDiscValues : NULL input mesh !");
812 return mesh->getCoordinatesAndOwner();
815 void MEDCouplingFieldDiscretizationOnNodes::computeMeshRestrictionFromTupleIds(const MEDCouplingMesh *mesh, const int *tupleIdsBg, const int *tupleIdsEnd,
816 DataArrayInt *&cellRestriction, DataArrayInt *&trueTupleRestriction) const throw(INTERP_KERNEL::Exception)
819 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationOnNodes::computeMeshRestrictionFromTupleIds : NULL input mesh !");
820 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret1=mesh->getCellIdsFullyIncludedInNodeIds(tupleIdsBg,tupleIdsEnd);
821 const MEDCouplingUMesh *meshc=dynamic_cast<const MEDCouplingUMesh *>(mesh);
823 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationOnNodes::computeMeshRestrictionFromTupleIds : trying to subpart field on nodes by node ids ! Your mesh has to be unstructured !");
824 MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> meshPart=static_cast<MEDCouplingUMesh *>(meshc->buildPartOfMySelf(ret1->begin(),ret1->end(),true));
825 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret2=meshPart->computeFetchedNodeIds();
826 cellRestriction=ret1.retn();
827 trueTupleRestriction=ret2.retn();
830 void MEDCouplingFieldDiscretizationOnNodes::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArray *da) const
833 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::checkCoherencyBetween : NULL input mesh or DataArray !");
834 if(mesh->getNumberOfNodes()!=da->getNumberOfTuples())
836 std::ostringstream message;
837 message << "Field on nodes invalid because there are " << mesh->getNumberOfNodes();
838 message << " nodes in mesh and " << da->getNumberOfTuples() << " tuples in field !";
839 throw INTERP_KERNEL::Exception(message.str().c_str());
844 * This method returns a submesh of 'mesh' instance constituting cell ids contained in array defined as an interval [start;end).
845 * @param di is an array returned that specifies entity ids (here nodes ids) in mesh 'mesh' of entity in returned submesh.
846 * Example : The first node id of returned mesh has the (*di)[0] id in 'mesh'
848 MEDCouplingMesh *MEDCouplingFieldDiscretizationOnNodes::buildSubMeshData(const MEDCouplingMesh *mesh, const int *start, const int *end, DataArrayInt *&di) const
851 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::buildSubMeshData : NULL input mesh !");
852 DataArrayInt *diTmp=0;
853 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPartAndReduceNodes(start,end,diTmp);
854 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> diTmpSafe(diTmp);
855 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> di2=diTmpSafe->invertArrayO2N2N2O(ret->getNumberOfNodes());
861 * This method is strictly equivalent to MEDCouplingFieldDiscretizationNodes::buildSubMeshData except that it is optimized for input defined as a range of cell ids.
863 * \param [out] beginOut Valid only if \a di is NULL
864 * \param [out] endOut Valid only if \a di is NULL
865 * \param [out] stepOut Valid only if \a di is NULL
866 * \param [out] di is an array returned that specifies entity ids (nodes, cells, Gauss points... ) in array if no output range is foundable.
868 * \sa MEDCouplingFieldDiscretizationNodes::buildSubMeshData
870 MEDCouplingMesh *MEDCouplingFieldDiscretizationOnNodes::buildSubMeshDataRange(const MEDCouplingMesh *mesh, int beginCellIds, int endCellIds, int stepCellIds, int& beginOut, int& endOut, int& stepOut, DataArrayInt *&di) const
873 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationOnNodes::buildSubMeshDataRange : NULL input mesh !");
874 DataArrayInt *diTmp=0;
875 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPartRangeAndReduceNodes(beginCellIds,endCellIds,stepCellIds,beginOut,endOut,stepOut,diTmp);
878 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> diTmpSafe(diTmp);
879 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> di2=diTmpSafe->invertArrayO2N2N2O(ret->getNumberOfNodes());
886 * This method returns a tuple ids selection from cell ids selection [start;end).
887 * This method is called by MEDCouplingFieldDiscretizationOnNodes::buildSubMeshData to return parameter \b di.
888 * Here for P1 only nodes fetched by submesh of mesh[startCellIds:endCellIds) is returned !
890 * \return a newly allocated array containing ids to select into the DataArrayDouble of the field.
893 DataArrayInt *MEDCouplingFieldDiscretizationOnNodes::computeTupleIdsToSelectFromCellIds(const MEDCouplingMesh *mesh, const int *startCellIds, const int *endCellIds) const
896 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP1::computeTupleIdsToSelectFromCellIds : NULL input mesh !");
897 const MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> umesh=mesh->buildUnstructured();
898 MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> umesh2=static_cast<MEDCouplingUMesh *>(umesh->buildPartOfMySelf(startCellIds,endCellIds,true));
899 return umesh2->computeFetchedNodeIds();
902 void MEDCouplingFieldDiscretizationOnNodes::renumberValuesOnNodes(double epsOnVals, const int *old2NewPtr, int newNbOfNodes, DataArrayDouble *arr) const
904 RenumberEntitiesFromO2NArr(epsOnVals,old2NewPtr,newNbOfNodes,arr,"Node");
908 * Nothing to do it's not a bug.
910 void MEDCouplingFieldDiscretizationOnNodes::renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const
915 * Nothing to do it's not a bug.
917 void MEDCouplingFieldDiscretizationOnNodes::renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const
921 void MEDCouplingFieldDiscretizationOnNodes::getValueOnPos(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, int i, int j, int k, double *res) const
923 const MEDCouplingCMesh *meshC=dynamic_cast<const MEDCouplingCMesh *>(mesh);
925 throw INTERP_KERNEL::Exception("OnNodes::getValueOnPos(i,j,k) is only accessible for structured meshes !");
926 int id=meshC->getNodeIdFromPos(i,j,k);
927 arr->getTuple(id,res);
930 TypeOfField MEDCouplingFieldDiscretizationP1::getEnum() const
936 * This method is simply called by MEDCouplingFieldDiscretization::deepCpy. It performs the deep copy of \a this.
938 * \sa MEDCouplingFieldDiscretization::deepCpy.
940 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationP1::clone() const
942 return new MEDCouplingFieldDiscretizationP1;
945 std::string MEDCouplingFieldDiscretizationP1::getStringRepr() const
947 return std::string(REPR);
950 const char *MEDCouplingFieldDiscretizationP1::getRepr() const
955 bool MEDCouplingFieldDiscretizationP1::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
959 reason="other spatial discretization is NULL, and this spatial discretization (P1) is defined.";
962 const MEDCouplingFieldDiscretizationP1 *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationP1 *>(other);
965 reason="Spatial discrtization of this is ON_NODES, which is not the case of other.";
969 void MEDCouplingFieldDiscretizationP1::checkCompatibilityWithNature(NatureOfField nat) const
971 if(nat!=ConservativeVolumic)
972 throw INTERP_KERNEL::Exception("Invalid nature for P1 field : expected ConservativeVolumic !");
975 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationP1::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
978 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP1::getMeasureField : mesh instance specified is NULL !");
979 return mesh->getMeasureFieldOnNode(isAbs);
982 void MEDCouplingFieldDiscretizationP1::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
985 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP1::getValueOn : NULL input mesh !");
986 int id=mesh->getCellContainingPoint(loc,_precision);
988 throw INTERP_KERNEL::Exception("Specified point is detected outside of mesh : unable to apply P1::getValueOn !");
989 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(id);
990 if(type!=INTERP_KERNEL::NORM_SEG2 && type!=INTERP_KERNEL::NORM_TRI3 && type!=INTERP_KERNEL::NORM_TETRA4)
991 throw INTERP_KERNEL::Exception("P1 getValueOn is not specified for not simplex cells !");
992 getValueInCell(mesh,id,arr,loc,res);
996 * This method localizes a point defined by 'loc' in a cell with id 'cellId' into mesh 'mesh'.
997 * The result is put into res expected to be of size at least arr->getNumberOfComponents()
999 void MEDCouplingFieldDiscretizationP1::getValueInCell(const MEDCouplingMesh *mesh, int cellId, const DataArrayDouble *arr, const double *loc, double *res) const
1002 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP1::getValueInCell : NULL input mesh !");
1003 std::vector<int> conn;
1004 std::vector<double> coo;
1005 mesh->getNodeIdsOfCell(cellId,conn);
1006 for(std::vector<int>::const_iterator iter=conn.begin();iter!=conn.end();iter++)
1007 mesh->getCoordinatesOfNode(*iter,coo);
1008 int spaceDim=mesh->getSpaceDimension();
1009 std::size_t nbOfNodes=conn.size();
1010 std::vector<const double *> vec(nbOfNodes);
1011 for(std::size_t i=0;i<nbOfNodes;i++)
1012 vec[i]=&coo[i*spaceDim];
1013 INTERP_KERNEL::AutoPtr<double> tmp=new double[nbOfNodes];
1014 INTERP_KERNEL::barycentric_coords(vec,loc,tmp);
1015 int sz=arr->getNumberOfComponents();
1016 INTERP_KERNEL::AutoPtr<double> tmp2=new double[sz];
1017 std::fill(res,res+sz,0.);
1018 for(std::size_t i=0;i<nbOfNodes;i++)
1020 arr->getTuple(conn[i],(double *)tmp2);
1021 std::transform((double *)tmp2,((double *)tmp2)+sz,(double *)tmp2,std::bind2nd(std::multiplies<double>(),tmp[i]));
1022 std::transform(res,res+sz,(double *)tmp2,res,std::plus<double>());
1026 DataArrayDouble *MEDCouplingFieldDiscretizationP1::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const
1029 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP1::getValueOnMulti : NULL input mesh !");
1030 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> eltsArr,eltsIndexArr;
1031 mesh->getCellsContainingPoints(loc,nbOfPoints,_precision,eltsArr,eltsIndexArr);
1032 const int *elts(eltsArr->begin()),*eltsIndex(eltsIndexArr->begin());
1033 int spaceDim=mesh->getSpaceDimension();
1034 int nbOfComponents=arr->getNumberOfComponents();
1035 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
1036 ret->alloc(nbOfPoints,nbOfComponents);
1037 double *ptToFill=ret->getPointer();
1038 for(int i=0;i<nbOfPoints;i++)
1039 if(eltsIndex[i+1]-eltsIndex[i]>=1)
1040 getValueInCell(mesh,elts[eltsIndex[i]],arr,loc+i*spaceDim,ptToFill+i*nbOfComponents);
1043 std::ostringstream oss; oss << "Point #" << i << " with coordinates : (";
1044 std::copy(loc+i*spaceDim,loc+(i+1)*spaceDim,std::ostream_iterator<double>(oss,", "));
1045 oss << ") detected outside mesh : unable to apply P1::getValueOnMulti ! ";
1046 throw INTERP_KERNEL::Exception(oss.str().c_str());
1051 void MEDCouplingFieldDiscretizationP1::reprQuickOverview(std::ostream& stream) const
1053 stream << "P1 spatial discretization.";
1056 MEDCouplingFieldDiscretizationPerCell::MEDCouplingFieldDiscretizationPerCell():_discr_per_cell(0)
1060 MEDCouplingFieldDiscretizationPerCell::~MEDCouplingFieldDiscretizationPerCell()
1063 _discr_per_cell->decrRef();
1067 * This constructor deep copies ParaMEDMEM::DataArrayInt instance from other (if any).
1069 MEDCouplingFieldDiscretizationPerCell::MEDCouplingFieldDiscretizationPerCell(const MEDCouplingFieldDiscretizationPerCell& other, const int *startCellIds, const int *endCellIds):_discr_per_cell(0)
1071 DataArrayInt *arr=other._discr_per_cell;
1074 if(startCellIds==0 && endCellIds==0)
1075 _discr_per_cell=arr->deepCpy();
1077 _discr_per_cell=arr->selectByTupleIdSafe(startCellIds,endCellIds);
1081 MEDCouplingFieldDiscretizationPerCell::MEDCouplingFieldDiscretizationPerCell(const MEDCouplingFieldDiscretizationPerCell& other, int beginCellIds, int endCellIds, int stepCellIds):_discr_per_cell(0)
1083 DataArrayInt *arr=other._discr_per_cell;
1086 _discr_per_cell=arr->selectByTupleId2(beginCellIds,endCellIds,stepCellIds);
1090 void MEDCouplingFieldDiscretizationPerCell::updateTime() const
1093 updateTimeWith(*_discr_per_cell);
1096 std::size_t MEDCouplingFieldDiscretizationPerCell::getHeapMemorySizeWithoutChildren() const
1098 std::size_t ret(MEDCouplingFieldDiscretization::getHeapMemorySizeWithoutChildren());
1102 std::vector<const BigMemoryObject *> MEDCouplingFieldDiscretizationPerCell::getDirectChildren() const
1104 std::vector<const BigMemoryObject *> ret(MEDCouplingFieldDiscretization::getDirectChildren());
1106 ret.push_back(_discr_per_cell);
1110 void MEDCouplingFieldDiscretizationPerCell::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArray *da) const
1112 if(!_discr_per_cell)
1113 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell has no discretization per cell !");
1115 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell::checkCoherencyBetween : NULL input mesh or DataArray !");
1116 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1117 if(nbOfTuples!=mesh->getNumberOfCells())
1118 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell has a discretization per cell but it's not matching the underlying mesh !");
1121 bool MEDCouplingFieldDiscretizationPerCell::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
1125 reason="other spatial discretization is NULL, and this spatial discretization (PerCell) is defined.";
1128 const MEDCouplingFieldDiscretizationPerCell *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationPerCell *>(other);
1131 reason="Spatial discretization of this is ON_GAUSS, which is not the case of other.";
1134 if(_discr_per_cell==0)
1135 return otherC->_discr_per_cell==0;
1136 if(otherC->_discr_per_cell==0)
1138 bool ret=_discr_per_cell->isEqualIfNotWhy(*otherC->_discr_per_cell,reason);
1140 reason.insert(0,"Field discretization per cell DataArrayInt given the discid per cell :");
1144 bool MEDCouplingFieldDiscretizationPerCell::isEqualWithoutConsideringStr(const MEDCouplingFieldDiscretization *other, double eps) const
1146 const MEDCouplingFieldDiscretizationPerCell *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationPerCell *>(other);
1149 if(_discr_per_cell==0)
1150 return otherC->_discr_per_cell==0;
1151 if(otherC->_discr_per_cell==0)
1153 return _discr_per_cell->isEqualWithoutConsideringStr(*otherC->_discr_per_cell);
1157 * This method is typically the first step of renumbering. The impact on _discr_per_cell is necessary here.
1158 * virtualy by this method.
1160 void MEDCouplingFieldDiscretizationPerCell::renumberCells(const int *old2NewBg, bool check)
1162 int nbCells=_discr_per_cell->getNumberOfTuples();
1163 const int *array=old2NewBg;
1165 array=DataArrayInt::CheckAndPreparePermutation(old2NewBg,old2NewBg+nbCells);
1167 DataArrayInt *dpc=_discr_per_cell->renumber(array);
1168 _discr_per_cell->decrRef();
1169 _discr_per_cell=dpc;
1172 free(const_cast<int *>(array));
1175 void MEDCouplingFieldDiscretizationPerCell::buildDiscrPerCellIfNecessary(const MEDCouplingMesh *mesh)
1178 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell::buildDiscrPerCellIfNecessary : NULL input mesh !");
1179 if(!_discr_per_cell)
1181 _discr_per_cell=DataArrayInt::New();
1182 int nbTuples=mesh->getNumberOfCells();
1183 _discr_per_cell->alloc(nbTuples,1);
1184 int *ptr=_discr_per_cell->getPointer();
1185 std::fill(ptr,ptr+nbTuples,DFT_INVALID_LOCID_VALUE);
1189 void MEDCouplingFieldDiscretizationPerCell::checkNoOrphanCells() const
1191 if(!_discr_per_cell)
1192 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell::checkNoOrphanCells : no discretization defined !");
1193 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> test=_discr_per_cell->getIdsEqual(DFT_INVALID_LOCID_VALUE);
1194 if(test->getNumberOfTuples()!=0)
1195 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell::checkNoOrphanCells : presence of orphan cells !");
1199 * This method is useful when 'this' describes a field discretization with several gauss discretization on a \b same cell type.
1200 * For example same NORM_TRI3 cells having 6 gauss points and others with 12 gauss points.
1201 * This method returns 2 arrays with same size : the return value and 'locIds' output parameter.
1202 * 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.
1203 * The return vector contains a set of newly created instance to deal with.
1204 * The returned vector represents a \b partition of cells ids with a gauss discretization set.
1206 * If no descretization is set in 'this' and exception will be thrown.
1208 std::vector<DataArrayInt *> MEDCouplingFieldDiscretizationPerCell::splitIntoSingleGaussDicrPerCellType(std::vector<int>& locIds) const
1210 if(!_discr_per_cell)
1211 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell::splitIntoSingleGaussDicrPerCellType : no descretization set !");
1212 return _discr_per_cell->partitionByDifferentValues(locIds);
1215 const DataArrayInt *MEDCouplingFieldDiscretizationPerCell::getArrayOfDiscIds() const
1217 return _discr_per_cell;
1220 void MEDCouplingFieldDiscretizationPerCell::setArrayOfDiscIds(const DataArrayInt *adids)
1222 if(adids!=_discr_per_cell)
1225 _discr_per_cell->decrRef();
1226 _discr_per_cell=const_cast<DataArrayInt *>(adids);
1228 _discr_per_cell->incrRef();
1233 MEDCouplingFieldDiscretizationGauss::MEDCouplingFieldDiscretizationGauss()
1237 MEDCouplingFieldDiscretizationGauss::MEDCouplingFieldDiscretizationGauss(const MEDCouplingFieldDiscretizationGauss& other, const int *startCellIds, const int *endCellIds):MEDCouplingFieldDiscretizationPerCell(other,startCellIds,endCellIds),_loc(other._loc)
1241 MEDCouplingFieldDiscretizationGauss::MEDCouplingFieldDiscretizationGauss(const MEDCouplingFieldDiscretizationGauss& other, int beginCellIds, int endCellIds, int stepCellIds):MEDCouplingFieldDiscretizationPerCell(other,beginCellIds,endCellIds,stepCellIds),_loc(other._loc)
1245 TypeOfField MEDCouplingFieldDiscretizationGauss::getEnum() const
1250 bool MEDCouplingFieldDiscretizationGauss::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
1254 reason="other spatial discretization is NULL, and this spatial discretization (Gauss) is defined.";
1257 const MEDCouplingFieldDiscretizationGauss *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationGauss *>(other);
1260 reason="Spatial discrtization of this is ON_GAUSS, which is not the case of other.";
1263 if(!MEDCouplingFieldDiscretizationPerCell::isEqualIfNotWhy(other,eps,reason))
1265 if(_loc.size()!=otherC->_loc.size())
1267 reason="Gauss spatial discretization : localization sizes differ";
1270 std::size_t sz=_loc.size();
1271 for(std::size_t i=0;i<sz;i++)
1272 if(!_loc[i].isEqual(otherC->_loc[i],eps))
1274 std::ostringstream oss; oss << "Gauss spatial discretization : Localization #" << i << " differ from this to other.";
1281 bool MEDCouplingFieldDiscretizationGauss::isEqualWithoutConsideringStr(const MEDCouplingFieldDiscretization *other, double eps) const
1283 const MEDCouplingFieldDiscretizationGauss *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationGauss *>(other);
1286 if(!MEDCouplingFieldDiscretizationPerCell::isEqualWithoutConsideringStr(other,eps))
1288 if(_loc.size()!=otherC->_loc.size())
1290 std::size_t sz=_loc.size();
1291 for(std::size_t i=0;i<sz;i++)
1292 if(!_loc[i].isEqual(otherC->_loc[i],eps))
1298 * This method is simply called by MEDCouplingFieldDiscretization::deepCpy. It performs the deep copy of \a this.
1300 * \sa MEDCouplingFieldDiscretization::deepCpy.
1302 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationGauss::clone() const
1304 return new MEDCouplingFieldDiscretizationGauss(*this);
1307 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationGauss::clonePart(const int *startCellIds, const int *endCellIds) const
1309 return new MEDCouplingFieldDiscretizationGauss(*this,startCellIds,endCellIds);
1312 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationGauss::clonePartRange(int beginCellIds, int endCellIds, int stepCellIds) const
1314 return new MEDCouplingFieldDiscretizationGauss(*this,beginCellIds,endCellIds,stepCellIds);
1317 std::string MEDCouplingFieldDiscretizationGauss::getStringRepr() const
1319 std::ostringstream oss; oss << REPR << "." << std::endl;
1322 if(_discr_per_cell->isAllocated())
1324 oss << "Discretization per cell : ";
1325 std::copy(_discr_per_cell->begin(),_discr_per_cell->end(),std::ostream_iterator<int>(oss,", "));
1329 oss << "Presence of " << _loc.size() << " localizations." << std::endl;
1331 for(std::vector<MEDCouplingGaussLocalization>::const_iterator it=_loc.begin();it!=_loc.end();it++,i++)
1333 oss << "+++++ Localization #" << i << " +++++" << std::endl;
1334 oss << (*it).getStringRepr();
1335 oss << "++++++++++" << std::endl;
1340 std::size_t MEDCouplingFieldDiscretizationGauss::getHeapMemorySizeWithoutChildren() const
1342 std::size_t ret(MEDCouplingFieldDiscretizationPerCell::getHeapMemorySizeWithoutChildren());
1343 ret+=_loc.capacity()*sizeof(MEDCouplingGaussLocalization);
1344 for(std::vector<MEDCouplingGaussLocalization>::const_iterator it=_loc.begin();it!=_loc.end();it++)
1345 ret+=(*it).getMemorySize();
1349 const char *MEDCouplingFieldDiscretizationGauss::getRepr() const
1355 * This method returns the number of tuples regarding exclusively the input code \b without \b using \b a \b mesh \b in \b input.
1356 * The input code coherency is also checked regarding spatial discretization of \a this.
1357 * If an incoherency is detected, an exception will be thrown. If the input code is coherent, the number of tuples expected is returned.
1358 * The number of tuples expected is equal to those to have a valid field lying on \a this and having a mesh fitting perfectly the input code (geometric type distribution).
1360 int MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode(const std::vector<int>& code, const std::vector<const DataArrayInt *>& idsPerType) const
1362 if(!_discr_per_cell || !_discr_per_cell->isAllocated() || _discr_per_cell->getNumberOfComponents()!=1)
1363 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode");
1364 if(code.size()%3!=0)
1365 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode : invalid input code !");
1366 int nbOfSplit=(int)idsPerType.size();
1367 int nbOfTypes=(int)code.size()/3;
1369 for(int i=0;i<nbOfTypes;i++)
1371 int nbOfEltInChunk=code[3*i+1];
1372 if(nbOfEltInChunk<0)
1373 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode : invalid input code ! presence of negative value in a type !");
1374 int pos=code[3*i+2];
1377 if(pos<0 || pos>=nbOfSplit)
1379 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode : input code points to pos " << pos << " in typeid " << i << " ! Should be in [0," << nbOfSplit << ") !";
1380 throw INTERP_KERNEL::Exception(oss.str().c_str());
1382 const DataArrayInt *ids(idsPerType[pos]);
1383 if(!ids || !ids->isAllocated() || ids->getNumberOfComponents()!=1 || ids->getNumberOfTuples()!=nbOfEltInChunk || ids->getMinValueInArray()<0)
1385 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode : input pfl chunck at pos " << pos << " should have " << i << " tuples and one component and with ids all >=0 !";
1386 throw INTERP_KERNEL::Exception(oss.str().c_str());
1389 ret+=nbOfEltInChunk;
1391 if(ret!=_discr_per_cell->getNumberOfTuples())
1393 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode : input code points to " << ret << " cells whereas discretization percell array lgth is " << _discr_per_cell->getNumberOfTuples() << " !";
1394 throw INTERP_KERNEL::Exception(oss.str().c_str());
1396 return getNumberOfTuples(0);//0 is not an error ! It is to be sure that input mesh is not used
1399 int MEDCouplingFieldDiscretizationGauss::getNumberOfTuples(const MEDCouplingMesh *) const
1402 if (_discr_per_cell == 0)
1403 throw INTERP_KERNEL::Exception("Discretization is not initialized!");
1404 const int *dcPtr=_discr_per_cell->getConstPointer();
1405 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1406 int maxSz=(int)_loc.size();
1407 for(const int *w=dcPtr;w!=dcPtr+nbOfTuples;w++)
1409 if(*w>=0 && *w<maxSz)
1410 ret+=_loc[*w].getNumberOfGaussPt();
1413 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getNumberOfTuples : At cell #" << std::distance(dcPtr,w) << " localization id is " << *w << " should be in [0," << maxSz << ") !";
1414 throw INTERP_KERNEL::Exception(oss.str().c_str());
1420 int MEDCouplingFieldDiscretizationGauss::getNumberOfMeshPlaces(const MEDCouplingMesh *mesh) const
1423 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getNumberOfMeshPlaces : NULL input mesh !");
1424 return mesh->getNumberOfCells();
1428 * This method is redevelopped for performance reasons, but it is equivalent to a call to MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField
1429 * and a call to DataArrayDouble::computeOffsets2 on the returned array.
1431 DataArrayInt *MEDCouplingFieldDiscretizationGauss::getOffsetArr(const MEDCouplingMesh *mesh) const
1434 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getOffsetArr : NULL input mesh !");
1435 int nbOfTuples=mesh->getNumberOfCells();
1436 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret=DataArrayInt::New();
1437 ret->alloc(nbOfTuples+1,1);
1438 int *retPtr=ret->getPointer();
1439 const int *start=_discr_per_cell->getConstPointer();
1440 if(_discr_per_cell->getNumberOfTuples()!=nbOfTuples)
1441 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getOffsetArr : mismatch between the mesh and the discretization ids array length !");
1442 int maxPossible=(int)_loc.size();
1444 for(int i=0;i<nbOfTuples;i++,start++)
1446 if(*start>=0 && *start<maxPossible)
1447 retPtr[i+1]=retPtr[i]+_loc[*start].getNumberOfGaussPt();
1450 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getOffsetArr : At position #" << i << " the locid = " << *start << " whereas it should be in [0," << maxPossible << ") !";
1451 throw INTERP_KERNEL::Exception(oss.str().c_str());
1457 void MEDCouplingFieldDiscretizationGauss::renumberArraysForCell(const MEDCouplingMesh *mesh, const std::vector<DataArray *>& arrays,
1458 const int *old2NewBg, bool check) throw(INTERP_KERNEL::Exception)
1461 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::renumberArraysForCell : NULL input mesh !");
1462 const int *array=old2NewBg;
1464 array=DataArrayInt::CheckAndPreparePermutation(old2NewBg,old2NewBg+mesh->getNumberOfCells());
1465 int nbOfCells=_discr_per_cell->getNumberOfTuples();
1466 int nbOfTuples=getNumberOfTuples(0);
1467 const int *dcPtr=_discr_per_cell->getConstPointer();
1468 int *array2=new int[nbOfTuples];//stores the final conversion array old2New to give to arrays in renumberInPlace.
1469 int *array3=new int[nbOfCells];//store for each cell in present dcp array (already renumbered) the offset needed by each cell in new numbering.
1471 for(int i=1;i<nbOfCells;i++)
1472 array3[i]=array3[i-1]+_loc[dcPtr[i-1]].getNumberOfGaussPt();
1474 for(int i=0;i<nbOfCells;i++)
1476 int nbOfGaussPt=_loc[dcPtr[array[i]]].getNumberOfGaussPt();
1477 for(int k=0;k<nbOfGaussPt;k++,j++)
1478 array2[j]=array3[array[i]]+k;
1481 for(std::vector<DataArray *>::const_iterator it=arrays.begin();it!=arrays.end();it++)
1483 (*it)->renumberInPlace(array2);
1486 free(const_cast<int*>(array));
1489 DataArrayDouble *MEDCouplingFieldDiscretizationGauss::getLocalizationOfDiscValues(const MEDCouplingMesh *mesh) const
1492 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getLocalizationOfDiscValues : NULL input mesh !");
1493 checkNoOrphanCells();
1494 MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> umesh=mesh->buildUnstructured();//in general do nothing
1495 int nbOfTuples=getNumberOfTuples(mesh);
1496 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
1497 int spaceDim=mesh->getSpaceDimension();
1498 ret->alloc(nbOfTuples,spaceDim);
1499 std::vector< int > locIds;
1500 std::vector<DataArrayInt *> parts=splitIntoSingleGaussDicrPerCellType(locIds);
1501 std::vector< MEDCouplingAutoRefCountObjectPtr<DataArrayInt> > parts2(parts.size());
1502 std::copy(parts.begin(),parts.end(),parts2.begin());
1503 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> offsets=buildNbOfGaussPointPerCellField();
1504 offsets->computeOffsets();
1505 const int *ptrOffsets=offsets->getConstPointer();
1506 const double *coords=umesh->getCoords()->getConstPointer();
1507 const int *connI=umesh->getNodalConnectivityIndex()->getConstPointer();
1508 const int *conn=umesh->getNodalConnectivity()->getConstPointer();
1509 double *valsToFill=ret->getPointer();
1510 for(std::size_t i=0;i<parts2.size();i++)
1512 INTERP_KERNEL::GaussCoords calculator;
1514 const MEDCouplingGaussLocalization& cli=_loc[locIds[i]];//curLocInfo
1515 INTERP_KERNEL::NormalizedCellType typ=cli.getType();
1516 const std::vector<double>& wg=cli.getWeights();
1517 calculator.addGaussInfo(typ,INTERP_KERNEL::CellModel::GetCellModel(typ).getDimension(),
1518 &cli.getGaussCoords()[0],(int)wg.size(),&cli.getRefCoords()[0],
1519 INTERP_KERNEL::CellModel::GetCellModel(typ).getNumberOfNodes());
1521 int nbt=parts2[i]->getNumberOfTuples();
1522 for(const int *w=parts2[i]->getConstPointer();w!=parts2[i]->getConstPointer()+nbt;w++)
1523 calculator.calculateCoords(cli.getType(),coords,spaceDim,conn+connI[*w]+1,valsToFill+spaceDim*(ptrOffsets[*w]));
1525 ret->copyStringInfoFrom(*umesh->getCoords());
1529 void MEDCouplingFieldDiscretizationGauss::computeMeshRestrictionFromTupleIds(const MEDCouplingMesh *mesh, const int *tupleIdsBg, const int *tupleIdsEnd,
1530 DataArrayInt *&cellRestriction, DataArrayInt *&trueTupleRestriction) const throw(INTERP_KERNEL::Exception)
1533 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::computeMeshRestrictionFromTupleIds : NULL input mesh !");
1534 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> tmp=DataArrayInt::New(); tmp->alloc((int)std::distance(tupleIdsBg,tupleIdsEnd),1);
1535 std::copy(tupleIdsBg,tupleIdsEnd,tmp->getPointer());
1537 tmp=tmp->buildUnique();
1538 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=buildNbOfGaussPointPerCellField();
1539 nbOfNodesPerCell->computeOffsets2();
1540 nbOfNodesPerCell->searchRangesInListOfIds(tmp,cellRestriction,trueTupleRestriction);
1546 void MEDCouplingFieldDiscretizationGauss::checkCompatibilityWithNature(NatureOfField nat) const
1550 void MEDCouplingFieldDiscretizationGauss::getTinySerializationIntInformation(std::vector<int>& tinyInfo) const
1554 val=_discr_per_cell->getNumberOfTuples();
1555 tinyInfo.push_back(val);
1556 tinyInfo.push_back((int)_loc.size());
1558 tinyInfo.push_back(-1);
1560 tinyInfo.push_back(_loc[0].getDimension());
1561 for(std::vector<MEDCouplingGaussLocalization>::const_iterator iter=_loc.begin();iter!=_loc.end();iter++)
1562 (*iter).pushTinySerializationIntInfo(tinyInfo);
1565 void MEDCouplingFieldDiscretizationGauss::getTinySerializationDbleInformation(std::vector<double>& tinyInfo) const
1567 for(std::vector<MEDCouplingGaussLocalization>::const_iterator iter=_loc.begin();iter!=_loc.end();iter++)
1568 (*iter).pushTinySerializationDblInfo(tinyInfo);
1571 void MEDCouplingFieldDiscretizationGauss::getSerializationIntArray(DataArrayInt *& arr) const
1575 arr=_discr_per_cell;
1578 void MEDCouplingFieldDiscretizationGauss::resizeForUnserialization(const std::vector<int>& tinyInfo, DataArrayInt *& arr)
1580 int val=tinyInfo[0];
1583 _discr_per_cell=DataArrayInt::New();
1584 _discr_per_cell->alloc(val,1);
1588 arr=_discr_per_cell;
1589 int nbOfLoc=tinyInfo[1];
1591 int dim=tinyInfo[2];
1594 delta=((int)tinyInfo.size()-3)/nbOfLoc;
1595 for(int i=0;i<nbOfLoc;i++)
1597 std::vector<int> tmp(tinyInfo.begin()+3+i*delta,tinyInfo.begin()+3+(i+1)*delta);
1598 MEDCouplingGaussLocalization elt=MEDCouplingGaussLocalization::BuildNewInstanceFromTinyInfo(dim,tmp);
1599 _loc.push_back(elt);
1603 void MEDCouplingFieldDiscretizationGauss::finishUnserialization(const std::vector<double>& tinyInfo)
1605 double *tmp=new double[tinyInfo.size()];
1606 std::copy(tinyInfo.begin(),tinyInfo.end(),tmp);
1607 const double *work=tmp;
1608 for(std::vector<MEDCouplingGaussLocalization>::iterator iter=_loc.begin();iter!=_loc.end();iter++)
1609 work=(*iter).fillWithValues(work);
1613 double MEDCouplingFieldDiscretizationGauss::getIJK(const MEDCouplingMesh *mesh, const DataArrayDouble *da, int cellId, int nodeIdInCell, int compoId) const
1615 int offset=getOffsetOfCell(cellId);
1616 return da->getIJ(offset+nodeIdInCell,compoId);
1619 void MEDCouplingFieldDiscretizationGauss::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArray *da) const
1622 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::checkCoherencyBetween : NULL input mesh or DataArray !");
1623 MEDCouplingFieldDiscretizationPerCell::checkCoherencyBetween(mesh,da);
1624 for(std::vector<MEDCouplingGaussLocalization>::const_iterator iter=_loc.begin();iter!=_loc.end();iter++)
1625 (*iter).checkCoherency();
1626 int nbOfDesc=(int)_loc.size();
1627 int nbOfCells=mesh->getNumberOfCells();
1628 const int *dc=_discr_per_cell->getConstPointer();
1629 for(int i=0;i<nbOfCells;i++)
1633 std::ostringstream oss; oss << "Cell # " << i << " of mesh \"" << mesh->getName() << "\" has an undefined gauss location ! Should never happend !";
1634 throw INTERP_KERNEL::Exception(oss.str().c_str());
1638 std::ostringstream oss; oss << "Cell # " << i << " of mesh \"" << mesh->getName() << "\" has no gauss location !";
1639 throw INTERP_KERNEL::Exception(oss.str().c_str());
1641 if(mesh->getTypeOfCell(i)!=_loc[dc[i]].getType())
1643 std::ostringstream oss; oss << "Types of mesh and gauss location mismatch for cell # " << i;
1644 throw INTERP_KERNEL::Exception(oss.str().c_str());
1647 int nbOfTuples=getNumberOfTuples(mesh);
1648 if(nbOfTuples!=da->getNumberOfTuples())
1650 std::ostringstream oss; oss << "Invalid number of tuples in the array : expecting " << nbOfTuples << " !";
1651 throw INTERP_KERNEL::Exception(oss.str().c_str());
1655 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationGauss::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
1658 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getMeasureField : mesh instance specified is NULL !");
1659 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=mesh->getMeasureField(isAbs);
1660 const double *volPtr=vol->getArray()->begin();
1661 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> ret=MEDCouplingFieldDouble::New(ON_GAUSS_PT);
1663 ret->setDiscretization(const_cast<MEDCouplingFieldDiscretizationGauss *>(this));
1664 if(!_discr_per_cell)
1665 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getMeasureField : no discr per cell array not defined ! spatial localization is incorrect !");
1666 _discr_per_cell->checkAllocated();
1667 if(_discr_per_cell->getNumberOfComponents()!=1)
1668 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getMeasureField : no discr per cell array defined but with nb of components different from 1 !");
1669 if(_discr_per_cell->getNumberOfTuples()!=vol->getNumberOfTuples())
1670 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 !");
1671 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> offset=getOffsetArr(mesh);
1672 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arr=DataArrayDouble::New(); arr->alloc(getNumberOfTuples(mesh),1);
1674 double *arrPtr=arr->getPointer();
1675 const int *offsetPtr=offset->getConstPointer();
1676 int maxGaussLoc=(int)_loc.size();
1677 std::vector<int> locIds;
1678 std::vector<DataArrayInt *> ids=splitIntoSingleGaussDicrPerCellType(locIds);
1679 std::vector< MEDCouplingAutoRefCountObjectPtr<DataArrayInt> > ids2(ids.size()); std::copy(ids.begin(),ids.end(),ids2.begin());
1680 for(std::size_t i=0;i<locIds.size();i++)
1682 const DataArrayInt *curIds=ids[i];
1683 int locId=locIds[i];
1684 if(locId>=0 && locId<maxGaussLoc)
1686 const MEDCouplingGaussLocalization& loc=_loc[locId];
1687 int nbOfGaussPt=loc.getNumberOfGaussPt();
1688 INTERP_KERNEL::AutoPtr<double> weights=new double[nbOfGaussPt];
1689 double sum=std::accumulate(loc.getWeights().begin(),loc.getWeights().end(),0.);
1690 std::transform(loc.getWeights().begin(),loc.getWeights().end(),(double *)weights,std::bind2nd(std::multiplies<double>(),1./sum));
1691 for(const int *cellId=curIds->begin();cellId!=curIds->end();cellId++)
1692 for(int j=0;j<nbOfGaussPt;j++)
1693 arrPtr[offsetPtr[*cellId]+j]=weights[j]*volPtr[*cellId];
1697 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getMeasureField : Presence of localization id " << locId << " in cell #" << curIds->getIJ(0,0) << " ! Must be in [0," << maxGaussLoc << ") !";
1698 throw INTERP_KERNEL::Exception(oss.str().c_str());
1701 ret->synchronizeTimeWithSupport();
1705 void MEDCouplingFieldDiscretizationGauss::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
1707 throw INTERP_KERNEL::Exception("Not implemented yet !");
1710 void MEDCouplingFieldDiscretizationGauss::getValueOnPos(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, int i, int j, int k, double *res) const
1712 throw INTERP_KERNEL::Exception("getValueOnPos(i,j,k) : Not applyable for Gauss points !");
1715 DataArrayDouble *MEDCouplingFieldDiscretizationGauss::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const
1717 throw INTERP_KERNEL::Exception("getValueOnMulti : Not implemented yet for gauss points !");
1720 MEDCouplingMesh *MEDCouplingFieldDiscretizationGauss::buildSubMeshData(const MEDCouplingMesh *mesh, const int *start, const int *end, DataArrayInt *&di) const
1723 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::buildSubMeshData : NULL input mesh !");
1724 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> diSafe=computeTupleIdsToSelectFromCellIds(mesh,start,end);
1725 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPart(start,end);
1731 * This method is strictly equivalent to MEDCouplingFieldDiscretizationGauss::buildSubMeshData except that it is optimized for input defined as a range of cell ids.
1733 * \param [out] beginOut Valid only if \a di is NULL
1734 * \param [out] endOut Valid only if \a di is NULL
1735 * \param [out] stepOut Valid only if \a di is NULL
1736 * \param [out] di is an array returned that specifies entity ids (nodes, cells, Gauss points... ) in array if no output range is foundable.
1738 * \sa MEDCouplingFieldDiscretizationGauss::buildSubMeshData
1740 MEDCouplingMesh *MEDCouplingFieldDiscretizationGauss::buildSubMeshDataRange(const MEDCouplingMesh *mesh, int beginCellIds, int endCellIds, int stepCellIds, int& beginOut, int& endOut, int& stepOut, DataArrayInt *&di) const
1742 if(stepCellIds!=1)//even for stepCellIds==-1 the output will not be a range
1743 return MEDCouplingFieldDiscretization::buildSubMeshDataRange(mesh,beginCellIds,endCellIds,stepCellIds,beginOut,endOut,stepOut,di);
1745 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::buildSubMeshDataRange : NULL input mesh !");
1746 if(!_discr_per_cell)
1747 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::buildSubMeshDataRange : no discretization array set !");
1748 di=0; beginOut=0; endOut=0; stepOut=stepCellIds;
1749 const char msg[]="MEDCouplingFieldDiscretizationGauss::buildSubMeshDataRange : cell #";
1750 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1751 const int *w=_discr_per_cell->begin();
1752 int nbMaxOfLocId=(int)_loc.size();
1753 for(int i=0;i<nbOfTuples;i++,w++)
1755 if(*w!=DFT_INVALID_LOCID_VALUE)
1757 if(*w>=0 && *w<nbMaxOfLocId)
1759 int delta=_loc[*w].getNumberOfGaussPt();
1767 { std::ostringstream oss; oss << msg << i << " has invalid id (" << *w << ") ! Should be in [0," << nbMaxOfLocId << ") !"; throw INTERP_KERNEL::Exception(oss.str().c_str()); }
1770 { std::ostringstream oss; oss << msg << i << " is detected as orphan !"; throw INTERP_KERNEL::Exception(oss.str().c_str()); }
1772 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPartRange(beginCellIds,endCellIds,stepCellIds);
1777 * This method returns a tuple ids selection from cell ids selection [start;end).
1778 * This method is called by MEDCouplingFieldDiscretizationGauss::buildSubMeshData to return parameter \b di.
1780 * \return a newly allocated array containing ids to select into the DataArrayDouble of the field.
1783 DataArrayInt *MEDCouplingFieldDiscretizationGauss::computeTupleIdsToSelectFromCellIds(const MEDCouplingMesh *mesh, const int *startCellIds, const int *endCellIds) const
1786 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::computeTupleIdsToSelectFromCellIds : null mesh !");
1787 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=buildNbOfGaussPointPerCellField();//check of _discr_per_cell not NULL pointer
1788 int nbOfCells=mesh->getNumberOfCells();
1789 if(_discr_per_cell->getNumberOfTuples()!=nbOfCells)
1790 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::computeTupleIdsToSelectFromCellIds : mismatch of nb of tuples of cell ids array and number of cells !");
1791 nbOfNodesPerCell->computeOffsets2();
1792 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> sel=DataArrayInt::New(); sel->useArray(startCellIds,false,CPP_DEALLOC,(int)std::distance(startCellIds,endCellIds),1);
1793 return sel->buildExplicitArrByRanges(nbOfNodesPerCell);
1797 * No implementation needed !
1799 void MEDCouplingFieldDiscretizationGauss::renumberValuesOnNodes(double , const int *, int newNbOfNodes, DataArrayDouble *) const
1803 void MEDCouplingFieldDiscretizationGauss::renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const
1805 throw INTERP_KERNEL::Exception("Not implemented yet !");
1808 void MEDCouplingFieldDiscretizationGauss::renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const
1810 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 !");
1813 void MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnType(const MEDCouplingMesh *mesh, INTERP_KERNEL::NormalizedCellType type, const std::vector<double>& refCoo,
1814 const std::vector<double>& gsCoo, const std::vector<double>& wg) throw(INTERP_KERNEL::Exception)
1817 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnType : NULL input mesh !");
1818 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
1819 if((int)cm.getDimension()!=mesh->getMeshDimension())
1821 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnType : mismatch of dimensions ! MeshDim==" << mesh->getMeshDimension();
1822 oss << " whereas Type '" << cm.getRepr() << "' has dimension " << cm.getDimension() << " !";
1823 throw INTERP_KERNEL::Exception(oss.str().c_str());
1825 buildDiscrPerCellIfNecessary(mesh);
1826 int id=(int)_loc.size();
1827 MEDCouplingGaussLocalization elt(type,refCoo,gsCoo,wg);
1828 _loc.push_back(elt);
1829 int *ptr=_discr_per_cell->getPointer();
1830 int nbCells=mesh->getNumberOfCells();
1831 for(int i=0;i<nbCells;i++)
1832 if(mesh->getTypeOfCell(i)==type)
1834 zipGaussLocalizations();
1837 void MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnCells(const MEDCouplingMesh *mesh, const int *begin, const int *end, const std::vector<double>& refCoo,
1838 const std::vector<double>& gsCoo, const std::vector<double>& wg) throw(INTERP_KERNEL::Exception)
1841 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnCells : NULL input mesh !");
1842 buildDiscrPerCellIfNecessary(mesh);
1843 if(std::distance(begin,end)<1)
1844 throw INTERP_KERNEL::Exception("Size of [begin,end) must be equal or greater than 1 !");
1845 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(*begin);
1846 MEDCouplingGaussLocalization elt(type,refCoo,gsCoo,wg);
1847 int id=(int)_loc.size();
1848 int *ptr=_discr_per_cell->getPointer();
1849 for(const int *w=begin+1;w!=end;w++)
1851 if(mesh->getTypeOfCell(*w)!=type)
1853 std::ostringstream oss; oss << "The cell with id " << *w << " has been detected to be incompatible in the [begin,end) array specified !";
1854 throw INTERP_KERNEL::Exception(oss.str().c_str());
1858 for(const int *w2=begin;w2!=end;w2++)
1861 _loc.push_back(elt);
1862 zipGaussLocalizations();
1865 void MEDCouplingFieldDiscretizationGauss::clearGaussLocalizations()
1869 _discr_per_cell->decrRef();
1875 void MEDCouplingFieldDiscretizationGauss::setGaussLocalization(int locId, const MEDCouplingGaussLocalization& loc)
1878 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::setGaussLocalization : localization id has to be >=0 !");
1879 int sz=(int)_loc.size();
1880 MEDCouplingGaussLocalization gLoc(INTERP_KERNEL::NORM_ERROR);
1882 _loc.resize(locId+1,gLoc);
1886 void MEDCouplingFieldDiscretizationGauss::resizeLocalizationVector(int newSz)
1889 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::resizeLocalizationVector : new size has to be >=0 !");
1890 MEDCouplingGaussLocalization gLoc(INTERP_KERNEL::NORM_ERROR);
1891 _loc.resize(newSz,gLoc);
1894 MEDCouplingGaussLocalization& MEDCouplingFieldDiscretizationGauss::getGaussLocalization(int locId)
1896 checkLocalizationId(locId);
1900 int MEDCouplingFieldDiscretizationGauss::getNbOfGaussLocalization() const
1902 return (int)_loc.size();
1905 int MEDCouplingFieldDiscretizationGauss::getGaussLocalizationIdOfOneCell(int cellId) const
1907 if(!_discr_per_cell)
1908 throw INTERP_KERNEL::Exception("No Gauss localization still set !");
1909 int locId=_discr_per_cell->begin()[cellId];
1911 throw INTERP_KERNEL::Exception("No Gauss localization set for the specified cell !");
1915 int MEDCouplingFieldDiscretizationGauss::getGaussLocalizationIdOfOneType(INTERP_KERNEL::NormalizedCellType type) const
1917 std::set<int> ret=getGaussLocalizationIdsOfOneType(type);
1919 throw INTERP_KERNEL::Exception("No gauss discretization found for the specified type !");
1921 throw INTERP_KERNEL::Exception("Several gauss discretizations have been found for the specified type !");
1922 return *ret.begin();
1925 std::set<int> MEDCouplingFieldDiscretizationGauss::getGaussLocalizationIdsOfOneType(INTERP_KERNEL::NormalizedCellType type) const
1927 if(!_discr_per_cell)
1928 throw INTERP_KERNEL::Exception("No Gauss localization still set !");
1931 for(std::vector<MEDCouplingGaussLocalization>::const_iterator iter=_loc.begin();iter!=_loc.end();iter++,id++)
1932 if((*iter).getType()==type)
1937 void MEDCouplingFieldDiscretizationGauss::getCellIdsHavingGaussLocalization(int locId, std::vector<int>& cellIds) const
1939 if(locId<0 || locId>=(int)_loc.size())
1940 throw INTERP_KERNEL::Exception("Invalid locId given : must be in range [0:getNbOfGaussLocalization()) !");
1941 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1942 const int *ptr=_discr_per_cell->getConstPointer();
1943 for(int i=0;i<nbOfTuples;i++)
1945 cellIds.push_back(i);
1948 const MEDCouplingGaussLocalization& MEDCouplingFieldDiscretizationGauss::getGaussLocalization(int locId) const
1950 checkLocalizationId(locId);
1954 void MEDCouplingFieldDiscretizationGauss::checkLocalizationId(int locId) const
1956 if(locId<0 || locId>=(int)_loc.size())
1957 throw INTERP_KERNEL::Exception("Invalid locId given : must be in range [0:getNbOfGaussLocalization()) !");
1960 int MEDCouplingFieldDiscretizationGauss::getOffsetOfCell(int cellId) const
1963 const int *start=_discr_per_cell->getConstPointer();
1964 for(const int *w=start;w!=start+cellId;w++)
1965 ret+=_loc[*w].getNumberOfGaussPt();
1970 * This method do the assumption that there is no orphan cell. If there is an exception is thrown.
1971 * This method makes the assumption too that '_discr_per_cell' is defined. If not an exception is thrown.
1972 * This method returns a newly created array with number of tuples equals to '_discr_per_cell->getNumberOfTuples' and number of components equal to 1.
1973 * The i_th tuple in returned array is the number of gauss point if the corresponding cell.
1975 DataArrayInt *MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField() const
1977 if(!_discr_per_cell)
1978 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField : no discretization array set !");
1979 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1980 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret=DataArrayInt::New();
1981 const int *w=_discr_per_cell->begin();
1982 ret->alloc(nbOfTuples,1);
1983 int *valsToFill=ret->getPointer();
1984 int nbMaxOfLocId=(int)_loc.size();
1985 for(int i=0;i<nbOfTuples;i++,w++)
1986 if(*w!=DFT_INVALID_LOCID_VALUE)
1988 if(*w>=0 && *w<nbMaxOfLocId)
1989 valsToFill[i]=_loc[*w].getNumberOfGaussPt();
1992 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField : cell #" << i << " has invalid id (" << *w << ") ! Should be in [0," << nbMaxOfLocId << ") !";
1993 throw INTERP_KERNEL::Exception(oss.str().c_str());
1998 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField : cell #" << i << " is detected as orphan !";
1999 throw INTERP_KERNEL::Exception(oss.str().c_str());
2004 void MEDCouplingFieldDiscretizationGauss::reprQuickOverview(std::ostream& stream) const
2006 stream << "Gauss points spatial discretization.";
2010 * This method makes the assumption that _discr_per_cell is set.
2011 * This method reduces as much as possible number size of _loc.
2012 * This method is useful when several set on same cells has been done and that some Gauss Localization are no more used.
2014 void MEDCouplingFieldDiscretizationGauss::zipGaussLocalizations()
2016 const int *start=_discr_per_cell->begin();
2017 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
2018 INTERP_KERNEL::AutoPtr<int> tmp=new int[_loc.size()];
2019 std::fill((int *)tmp,(int *)tmp+_loc.size(),-2);
2020 for(const int *w=start;w!=start+nbOfTuples;w++)
2024 for(int i=0;i<(int)_loc.size();i++)
2027 if(fid==(int)_loc.size())
2030 int *start2=_discr_per_cell->getPointer();
2031 for(int *w2=start2;w2!=start2+nbOfTuples;w2++)
2034 std::vector<MEDCouplingGaussLocalization> tmpLoc;
2035 for(int i=0;i<(int)_loc.size();i++)
2037 tmpLoc.push_back(_loc[i]);
2041 MEDCouplingFieldDiscretizationGaussNE::MEDCouplingFieldDiscretizationGaussNE()
2045 TypeOfField MEDCouplingFieldDiscretizationGaussNE::getEnum() const
2051 * This method is simply called by MEDCouplingFieldDiscretization::deepCpy. It performs the deep copy of \a this.
2053 * \sa MEDCouplingFieldDiscretization::deepCpy.
2055 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationGaussNE::clone() const
2057 return new MEDCouplingFieldDiscretizationGaussNE(*this);
2060 std::string MEDCouplingFieldDiscretizationGaussNE::getStringRepr() const
2062 return std::string(REPR);
2065 const char *MEDCouplingFieldDiscretizationGaussNE::getRepr() const
2070 bool MEDCouplingFieldDiscretizationGaussNE::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
2074 reason="other spatial discretization is NULL, and this spatial discretization (GaussNE) is defined.";
2077 const MEDCouplingFieldDiscretizationGaussNE *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationGaussNE *>(other);
2080 reason="Spatial discrtization of this is ON_GAUSS_NE, which is not the case of other.";
2085 * This method returns the number of tuples regarding exclusively the input code \b without \b using \b a \b mesh \b in \b input.
2086 * The input code coherency is also checked regarding spatial discretization of \a this.
2087 * If an incoherency is detected, an exception will be thrown. If the input code is coherent, the number of tuples expected is returned.
2088 * The number of tuples expected is equal to those to have a valid field lying on \a this and having a mesh fitting perfectly the input code (geometric type distribution).
2090 int MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode(const std::vector<int>& code, const std::vector<const DataArrayInt *>& idsPerType) const
2092 if(code.size()%3!=0)
2093 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode : invalid input code !");
2094 int nbOfSplit=(int)idsPerType.size();
2095 int nbOfTypes=(int)code.size()/3;
2097 for(int i=0;i<nbOfTypes;i++)
2099 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel((INTERP_KERNEL::NormalizedCellType)code[3*i]));
2102 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode : At pos #" << i << " the geometric type " << cm.getRepr() << " is dynamic ! There are not managed by GAUSS_NE field discretization !";
2103 throw INTERP_KERNEL::Exception(oss.str().c_str());
2105 int nbOfEltInChunk=code[3*i+1];
2106 if(nbOfEltInChunk<0)
2107 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode : invalid input code ! presence of negative value in a type !");
2108 int pos=code[3*i+2];
2111 if(pos<0 || pos>=nbOfSplit)
2113 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode : input code points to pos " << pos << " in typeid " << i << " ! Should be in [0," << nbOfSplit << ") !";
2114 throw INTERP_KERNEL::Exception(oss.str().c_str());
2116 const DataArrayInt *ids(idsPerType[pos]);
2117 if(!ids || !ids->isAllocated() || ids->getNumberOfComponents()!=1 || ids->getNumberOfTuples()!=nbOfEltInChunk || ids->getMinValueInArray()<0)
2119 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode : input pfl chunck at pos " << pos << " should have " << i << " tuples and one component and with ids all >=0 !";
2120 throw INTERP_KERNEL::Exception(oss.str().c_str());
2123 ret+=nbOfEltInChunk*(int)cm.getNumberOfNodes();
2128 int MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuples(const MEDCouplingMesh *mesh) const
2131 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuples : NULL input mesh !");
2133 int nbOfCells=mesh->getNumberOfCells();
2134 for(int i=0;i<nbOfCells;i++)
2136 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
2137 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2139 throw INTERP_KERNEL::Exception("Not implemented yet Gauss node on elements for polygons and polyedrons !");
2140 ret+=cm.getNumberOfNodes();
2145 int MEDCouplingFieldDiscretizationGaussNE::getNumberOfMeshPlaces(const MEDCouplingMesh *mesh) const
2148 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getNumberOfMeshPlaces : NULL input mesh !");
2149 return mesh->getNumberOfCells();
2152 DataArrayInt *MEDCouplingFieldDiscretizationGaussNE::getOffsetArr(const MEDCouplingMesh *mesh) const
2155 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getOffsetArr : NULL input mesh !");
2156 int nbOfTuples=mesh->getNumberOfCells();
2157 DataArrayInt *ret=DataArrayInt::New();
2158 ret->alloc(nbOfTuples+1,1);
2159 int *retPtr=ret->getPointer();
2161 for(int i=0;i<nbOfTuples;i++)
2163 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
2164 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2166 throw INTERP_KERNEL::Exception("Not implemented yet Gauss node on elements for polygons and polyedrons !");
2167 retPtr[i+1]=retPtr[i]+cm.getNumberOfNodes();
2172 void MEDCouplingFieldDiscretizationGaussNE::renumberArraysForCell(const MEDCouplingMesh *mesh, const std::vector<DataArray *>& arrays,
2173 const int *old2NewBg, bool check) throw(INTERP_KERNEL::Exception)
2176 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::renumberArraysForCell : NULL input mesh !");
2177 const int *array=old2NewBg;
2179 array=DataArrayInt::CheckAndPreparePermutation(old2NewBg,old2NewBg+mesh->getNumberOfCells());
2180 int nbOfCells=mesh->getNumberOfCells();
2181 int nbOfTuples=getNumberOfTuples(mesh);
2182 int *array2=new int[nbOfTuples];//stores the final conversion array old2New to give to arrays in renumberInPlace.
2183 int *array3=new int[nbOfCells];//store for each cell in after renumbering the offset needed by each cell in new numbering.
2185 for(int i=1;i<nbOfCells;i++)
2187 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell((int)std::distance(array,std::find(array,array+nbOfCells,i-1)));
2188 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2189 array3[i]=array3[i-1]+cm.getNumberOfNodes();
2192 for(int i=0;i<nbOfCells;i++)
2194 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
2195 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2196 for(int k=0;k<(int)cm.getNumberOfNodes();k++,j++)
2197 array2[j]=array3[array[i]]+k;
2200 for(std::vector<DataArray *>::const_iterator it=arrays.begin();it!=arrays.end();it++)
2202 (*it)->renumberInPlace(array2);
2205 free(const_cast<int *>(array));
2208 DataArrayDouble *MEDCouplingFieldDiscretizationGaussNE::getLocalizationOfDiscValues(const MEDCouplingMesh *mesh) const
2211 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getLocalizationOfDiscValues : NULL input mesh !");
2212 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
2213 MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> umesh=mesh->buildUnstructured();//in general do nothing
2214 int nbOfTuples=getNumberOfTuples(umesh);
2215 int spaceDim=mesh->getSpaceDimension();
2216 ret->alloc(nbOfTuples,spaceDim);
2217 const double *coords=umesh->getCoords()->begin();
2218 const int *connI=umesh->getNodalConnectivityIndex()->getConstPointer();
2219 const int *conn=umesh->getNodalConnectivity()->getConstPointer();
2220 int nbCells=umesh->getNumberOfCells();
2221 double *retPtr=ret->getPointer();
2222 for(int i=0;i<nbCells;i++,connI++)
2223 for(const int *w=conn+connI[0]+1;w!=conn+connI[1];w++)
2225 retPtr=std::copy(coords+(*w)*spaceDim,coords+((*w)+1)*spaceDim,retPtr);
2230 * Reimplemented from MEDCouplingFieldDiscretization::integral for performance reason. The default implementation is valid too for GAUSS_NE spatial discretization.
2232 void MEDCouplingFieldDiscretizationGaussNE::integral(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, bool isWAbs, double *res) const
2235 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::integral : input mesh or array is null !");
2236 int nbOfCompo=arr->getNumberOfComponents();
2237 std::fill(res,res+nbOfCompo,0.);
2239 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=mesh->getMeasureField(isWAbs);
2240 std::set<INTERP_KERNEL::NormalizedCellType> types=mesh->getAllGeoTypes();
2241 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=mesh->computeNbOfNodesPerCell();
2242 nbOfNodesPerCell->computeOffsets2();
2243 const double *arrPtr=arr->begin(),*volPtr=vol->getArray()->begin();
2244 for(std::set<INTERP_KERNEL::NormalizedCellType>::const_iterator it=types.begin();it!=types.end();it++)
2246 std::size_t wArrSz=-1;
2247 const double *wArr=GetWeightArrayFromGeometricType(*it,wArrSz);
2248 INTERP_KERNEL::AutoPtr<double> wArr2=new double[wArrSz];
2249 double sum=std::accumulate(wArr,wArr+wArrSz,0.);
2250 std::transform(wArr,wArr+wArrSz,(double *)wArr2,std::bind2nd(std::multiplies<double>(),1./sum));
2251 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ids=mesh->giveCellsWithType(*it);
2252 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ids2=ids->buildExplicitArrByRanges(nbOfNodesPerCell);
2253 const int *ptIds2=ids2->begin(),*ptIds=ids->begin();
2254 int nbOfCellsWithCurGeoType=ids->getNumberOfTuples();
2255 for(int i=0;i<nbOfCellsWithCurGeoType;i++,ptIds++,ptIds2+=wArrSz)
2257 for(int k=0;k<nbOfCompo;k++)
2260 for(std::size_t j=0;j<wArrSz;j++)
2261 tmp+=arrPtr[nbOfCompo*ptIds2[j]+k]*wArr2[j];
2262 res[k]+=tmp*volPtr[*ptIds];
2268 const double *MEDCouplingFieldDiscretizationGaussNE::GetWeightArrayFromGeometricType(INTERP_KERNEL::NormalizedCellType geoType, std::size_t& lgth)
2272 case INTERP_KERNEL::NORM_SEG2:
2273 lgth=(int)sizeof(FGP_SEG2)/sizeof(double);
2275 case INTERP_KERNEL::NORM_SEG3:
2276 lgth=(int)sizeof(FGP_SEG3)/sizeof(double);
2278 case INTERP_KERNEL::NORM_SEG4:
2279 lgth=(int)sizeof(FGP_SEG4)/sizeof(double);
2281 case INTERP_KERNEL::NORM_TRI3:
2282 lgth=(int)sizeof(FGP_TRI3)/sizeof(double);
2284 case INTERP_KERNEL::NORM_TRI6:
2285 lgth=(int)sizeof(FGP_TRI6)/sizeof(double);
2287 case INTERP_KERNEL::NORM_TRI7:
2288 lgth=(int)sizeof(FGP_TRI7)/sizeof(double);
2290 case INTERP_KERNEL::NORM_QUAD4:
2291 lgth=(int)sizeof(FGP_QUAD4)/sizeof(double);
2293 case INTERP_KERNEL::NORM_QUAD8:
2294 lgth=(int)sizeof(FGP_QUAD8)/sizeof(double);
2296 case INTERP_KERNEL::NORM_QUAD9:
2297 lgth=(int)sizeof(FGP_QUAD9)/sizeof(double);
2299 case INTERP_KERNEL::NORM_TETRA4:
2300 lgth=(int)sizeof(FGP_TETRA4)/sizeof(double);
2302 case INTERP_KERNEL::NORM_PENTA6:
2303 lgth=(int)sizeof(FGP_PENTA6)/sizeof(double);
2305 case INTERP_KERNEL::NORM_HEXA8:
2306 lgth=(int)sizeof(FGP_HEXA8)/sizeof(double);
2308 case INTERP_KERNEL::NORM_HEXA27:
2309 lgth=(int)sizeof(FGP_HEXA27)/sizeof(double);
2311 case INTERP_KERNEL::NORM_PYRA5:
2312 lgth=(int)sizeof(FGP_PYRA5)/sizeof(double);
2315 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 !");
2319 const double *MEDCouplingFieldDiscretizationGaussNE::GetRefCoordsFromGeometricType(INTERP_KERNEL::NormalizedCellType geoType, std::size_t& lgth)
2323 case INTERP_KERNEL::NORM_SEG2:
2324 lgth=(int)sizeof(REF_SEG2)/sizeof(double);
2326 case INTERP_KERNEL::NORM_SEG3:
2327 lgth=(int)sizeof(REF_SEG3)/sizeof(double);
2329 case INTERP_KERNEL::NORM_SEG4:
2330 lgth=(int)sizeof(REF_SEG4)/sizeof(double);
2332 case INTERP_KERNEL::NORM_TRI3:
2333 lgth=(int)sizeof(REF_TRI3)/sizeof(double);
2335 case INTERP_KERNEL::NORM_TRI6:
2336 lgth=(int)sizeof(REF_TRI6)/sizeof(double);
2338 case INTERP_KERNEL::NORM_TRI7:
2339 lgth=(int)sizeof(REF_TRI7)/sizeof(double);
2341 case INTERP_KERNEL::NORM_QUAD4:
2342 lgth=(int)sizeof(REF_QUAD4)/sizeof(double);
2344 case INTERP_KERNEL::NORM_QUAD8:
2345 lgth=(int)sizeof(REF_QUAD8)/sizeof(double);
2347 case INTERP_KERNEL::NORM_QUAD9:
2348 lgth=(int)sizeof(REF_QUAD9)/sizeof(double);
2350 case INTERP_KERNEL::NORM_TETRA4:
2351 lgth=(int)sizeof(REF_TETRA4)/sizeof(double);
2353 case INTERP_KERNEL::NORM_TETRA10:
2354 lgth=(int)sizeof(REF_TETRA10)/sizeof(double);
2356 case INTERP_KERNEL::NORM_PENTA6:
2357 lgth=(int)sizeof(REF_PENTA6)/sizeof(double);
2359 case INTERP_KERNEL::NORM_PENTA15:
2360 lgth=(int)sizeof(REF_PENTA15)/sizeof(double);
2362 case INTERP_KERNEL::NORM_HEXA8:
2363 lgth=(int)sizeof(REF_HEXA8)/sizeof(double);
2365 case INTERP_KERNEL::NORM_HEXA20:
2366 lgth=(int)sizeof(REF_HEXA20)/sizeof(double);
2368 case INTERP_KERNEL::NORM_HEXA27:
2369 lgth=(int)sizeof(REF_HEXA27)/sizeof(double);
2371 case INTERP_KERNEL::NORM_PYRA5:
2372 lgth=(int)sizeof(REF_PYRA5)/sizeof(double);
2374 case INTERP_KERNEL::NORM_PYRA13:
2375 lgth=(int)sizeof(REF_PYRA13)/sizeof(double);
2378 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 !");
2382 const double *MEDCouplingFieldDiscretizationGaussNE::GetLocsFromGeometricType(INTERP_KERNEL::NormalizedCellType geoType, std::size_t& lgth)
2386 case INTERP_KERNEL::NORM_SEG2:
2388 lgth=(int)sizeof(LOC_SEG2)/sizeof(double);
2391 case INTERP_KERNEL::NORM_SEG3:
2393 lgth=(int)sizeof(LOC_SEG3)/sizeof(double);
2396 case INTERP_KERNEL::NORM_SEG4:
2398 lgth=(int)sizeof(LOC_SEG4)/sizeof(double);
2401 case INTERP_KERNEL::NORM_TRI3:
2403 lgth=(int)sizeof(LOC_TRI3)/sizeof(double);
2406 case INTERP_KERNEL::NORM_TRI6:
2408 lgth=(int)sizeof(LOC_TRI6)/sizeof(double);
2411 case INTERP_KERNEL::NORM_TRI7:
2413 lgth=(int)sizeof(LOC_TRI7)/sizeof(double);
2416 case INTERP_KERNEL::NORM_QUAD4:
2418 lgth=(int)sizeof(LOC_QUAD4)/sizeof(double);
2421 case INTERP_KERNEL::NORM_QUAD8:
2423 lgth=(int)sizeof(LOC_QUAD8)/sizeof(double);
2426 case INTERP_KERNEL::NORM_QUAD9:
2428 lgth=(int)sizeof(LOC_QUAD9)/sizeof(double);
2431 case INTERP_KERNEL::NORM_TETRA4:
2433 lgth=(int)sizeof(LOC_TETRA4)/sizeof(double);
2436 case INTERP_KERNEL::NORM_PENTA6:
2438 lgth=(int)sizeof(LOC_PENTA6)/sizeof(double);
2441 case INTERP_KERNEL::NORM_HEXA8:
2443 lgth=(int)sizeof(LOC_HEXA8)/sizeof(double);
2446 case INTERP_KERNEL::NORM_HEXA27:
2448 lgth=(int)sizeof(LOC_HEXA27)/sizeof(double);
2451 case INTERP_KERNEL::NORM_PYRA5:
2453 lgth=(int)sizeof(LOC_PYRA5)/sizeof(double);
2457 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::GetLocsFromGeometricType : 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 !");
2461 void MEDCouplingFieldDiscretizationGaussNE::computeMeshRestrictionFromTupleIds(const MEDCouplingMesh *mesh, const int *tupleIdsBg, const int *tupleIdsEnd,
2462 DataArrayInt *&cellRestriction, DataArrayInt *&trueTupleRestriction) const throw(INTERP_KERNEL::Exception)
2465 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::computeMeshRestrictionFromTupleIds : NULL input mesh !");
2466 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> tmp=DataArrayInt::New(); tmp->alloc((int)std::distance(tupleIdsBg,tupleIdsEnd),1);
2467 std::copy(tupleIdsBg,tupleIdsEnd,tmp->getPointer());
2469 tmp=tmp->buildUnique();
2470 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=mesh->computeNbOfNodesPerCell();
2471 nbOfNodesPerCell->computeOffsets2();
2472 nbOfNodesPerCell->searchRangesInListOfIds(tmp,cellRestriction,trueTupleRestriction);
2475 void MEDCouplingFieldDiscretizationGaussNE::checkCompatibilityWithNature(NatureOfField nat) const
2479 double MEDCouplingFieldDiscretizationGaussNE::getIJK(const MEDCouplingMesh *mesh, const DataArrayDouble *da, int cellId, int nodeIdInCell, int compoId) const
2482 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getIJK : NULL input mesh !");
2484 for(int i=0;i<cellId;i++)
2486 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
2487 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2488 offset+=cm.getNumberOfNodes();
2490 return da->getIJ(offset+nodeIdInCell,compoId);
2493 void MEDCouplingFieldDiscretizationGaussNE::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArray *da) const
2495 int nbOfTuples=getNumberOfTuples(mesh);
2496 if(nbOfTuples!=da->getNumberOfTuples())
2498 std::ostringstream oss; oss << "Invalid number of tuples in the array : expecting " << nbOfTuples << " !";
2499 throw INTERP_KERNEL::Exception(oss.str().c_str());
2503 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationGaussNE::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
2506 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getMeasureField : mesh instance specified is NULL !");
2507 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=mesh->getMeasureField(isAbs);
2508 const double *volPtr=vol->getArray()->begin();
2509 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> ret=MEDCouplingFieldDouble::New(ON_GAUSS_NE);
2512 std::set<INTERP_KERNEL::NormalizedCellType> types=mesh->getAllGeoTypes();
2513 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=mesh->computeNbOfNodesPerCell();
2514 int nbTuples=nbOfNodesPerCell->accumulate(0);
2515 nbOfNodesPerCell->computeOffsets2();
2516 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arr=DataArrayDouble::New(); arr->alloc(nbTuples,1);
2518 double *arrPtr=arr->getPointer();
2519 for(std::set<INTERP_KERNEL::NormalizedCellType>::const_iterator it=types.begin();it!=types.end();it++)
2521 std::size_t wArrSz=-1;
2522 const double *wArr=GetWeightArrayFromGeometricType(*it,wArrSz);
2523 INTERP_KERNEL::AutoPtr<double> wArr2=new double[wArrSz];
2524 double sum=std::accumulate(wArr,wArr+wArrSz,0.);
2525 std::transform(wArr,wArr+wArrSz,(double *)wArr2,std::bind2nd(std::multiplies<double>(),1./sum));
2526 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ids=mesh->giveCellsWithType(*it);
2527 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ids2=ids->buildExplicitArrByRanges(nbOfNodesPerCell);
2528 const int *ptIds2=ids2->begin(),*ptIds=ids->begin();
2529 int nbOfCellsWithCurGeoType=ids->getNumberOfTuples();
2530 for(int i=0;i<nbOfCellsWithCurGeoType;i++,ptIds++)
2531 for(std::size_t j=0;j<wArrSz;j++,ptIds2++)
2532 arrPtr[*ptIds2]=wArr2[j]*volPtr[*ptIds];
2534 ret->synchronizeTimeWithSupport();
2538 void MEDCouplingFieldDiscretizationGaussNE::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
2540 throw INTERP_KERNEL::Exception("Not implemented yet !");
2543 void MEDCouplingFieldDiscretizationGaussNE::getValueOnPos(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, int i, int j, int k, double *res) const
2545 throw INTERP_KERNEL::Exception("getValueOnPos(i,j,k) : Not applyable for Gauss points !");
2548 DataArrayDouble *MEDCouplingFieldDiscretizationGaussNE::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const
2550 throw INTERP_KERNEL::Exception("getValueOnMulti : Not implemented for Gauss NE !");
2553 MEDCouplingMesh *MEDCouplingFieldDiscretizationGaussNE::buildSubMeshData(const MEDCouplingMesh *mesh, const int *start, const int *end, DataArrayInt *&di) const
2556 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::buildSubMeshData : NULL input mesh !");
2557 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> diSafe=computeTupleIdsToSelectFromCellIds(mesh,start,end);
2558 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPart(start,end);
2564 * This method is strictly equivalent to MEDCouplingFieldDiscretizationGauss::buildSubMeshData except that it is optimized for input defined as a range of cell ids.
2566 * \param [out] beginOut Valid only if \a di is NULL
2567 * \param [out] endOut Valid only if \a di is NULL
2568 * \param [out] stepOut Valid only if \a di is NULL
2569 * \param [out] di is an array returned that specifies entity ids (nodes, cells, Gauss points... ) in array if no output range is foundable.
2571 * \sa MEDCouplingFieldDiscretizationGauss::buildSubMeshData
2573 MEDCouplingMesh *MEDCouplingFieldDiscretizationGaussNE::buildSubMeshDataRange(const MEDCouplingMesh *mesh, int beginCellIds, int endCellIds, int stepCellIds, int& beginOut, int& endOut, int& stepOut, DataArrayInt *&di) const
2575 if(stepCellIds!=1)//even for stepCellIds==-1 the output will not be a range
2576 return MEDCouplingFieldDiscretization::buildSubMeshDataRange(mesh,beginCellIds,endCellIds,stepCellIds,beginOut,endOut,stepOut,di);
2578 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::buildSubMeshDataRange : NULL input mesh !");
2579 int nbOfCells=mesh->getNumberOfCells();
2580 di=0; beginOut=0; endOut=0; stepOut=stepCellIds;
2581 const char msg[]="MEDCouplingFieldDiscretizationGaussNE::buildSubMeshDataRange : cell #";
2582 for(int i=0;i<nbOfCells;i++)
2584 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
2585 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2587 { std::ostringstream oss; oss << msg << i << " presence of dynamic cell (polygons and polyedrons) ! Not implemented !"; throw INTERP_KERNEL::Exception(oss.str().c_str()); }
2588 int delta=cm.getNumberOfNodes();
2595 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPartRange(beginCellIds,endCellIds,stepCellIds);
2601 * This method returns a tuple ids selection from cell ids selection [start;end).
2602 * This method is called by MEDCouplingFieldDiscretizationGaussNE::buildSubMeshData to return parameter \b di.
2604 * \return a newly allocated array containing ids to select into the DataArrayDouble of the field.
2607 DataArrayInt *MEDCouplingFieldDiscretizationGaussNE::computeTupleIdsToSelectFromCellIds(const MEDCouplingMesh *mesh, const int *startCellIds, const int *endCellIds) const
2610 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::computeTupleIdsToSelectFromCellIds : null mesh !");
2611 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=mesh->computeNbOfNodesPerCell();
2612 nbOfNodesPerCell->computeOffsets2();
2613 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> sel=DataArrayInt::New(); sel->useArray(startCellIds,false,CPP_DEALLOC,(int)std::distance(startCellIds,endCellIds),1);
2614 return sel->buildExplicitArrByRanges(nbOfNodesPerCell);
2618 * No implementation needed !
2620 void MEDCouplingFieldDiscretizationGaussNE::renumberValuesOnNodes(double , const int *, int newNbOfNodes, DataArrayDouble *) const
2624 void MEDCouplingFieldDiscretizationGaussNE::renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const
2626 throw INTERP_KERNEL::Exception("Not implemented yet !");
2629 void MEDCouplingFieldDiscretizationGaussNE::renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const
2631 throw INTERP_KERNEL::Exception("Not implemented yet !");
2634 void MEDCouplingFieldDiscretizationGaussNE::reprQuickOverview(std::ostream& stream) const
2636 stream << "Gauss points on nodes per element spatial discretization.";
2639 MEDCouplingFieldDiscretizationGaussNE::MEDCouplingFieldDiscretizationGaussNE(const MEDCouplingFieldDiscretizationGaussNE& other):MEDCouplingFieldDiscretization(other)
2643 TypeOfField MEDCouplingFieldDiscretizationKriging::getEnum() const
2648 const char *MEDCouplingFieldDiscretizationKriging::getRepr() const
2654 * This method is simply called by MEDCouplingFieldDiscretization::deepCpy. It performs the deep copy of \a this.
2656 * \sa MEDCouplingFieldDiscretization::deepCpy.
2658 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationKriging::clone() const
2660 return new MEDCouplingFieldDiscretizationKriging;
2663 std::string MEDCouplingFieldDiscretizationKriging::getStringRepr() const
2665 return std::string(REPR);
2668 void MEDCouplingFieldDiscretizationKriging::checkCompatibilityWithNature(NatureOfField nat) const
2670 if(nat!=ConservativeVolumic)
2671 throw INTERP_KERNEL::Exception("Invalid nature for Kriging field : expected ConservativeVolumic !");
2674 bool MEDCouplingFieldDiscretizationKriging::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
2678 reason="other spatial discretization is NULL, and this spatial discretization (Kriginig) is defined.";
2681 const MEDCouplingFieldDiscretizationKriging *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationKriging *>(other);
2684 reason="Spatial discrtization of this is ON_NODES_KR, which is not the case of other.";
2688 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationKriging::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
2691 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationKriging::getMeasureField : mesh instance specified is NULL !");
2692 throw INTERP_KERNEL::Exception("getMeasureField on FieldDiscretizationKriging : not implemented yet !");
2695 void MEDCouplingFieldDiscretizationKriging::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
2697 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> res2=MEDCouplingFieldDiscretizationKriging::getValueOnMulti(arr,mesh,loc,1);
2698 std::copy(res2->begin(),res2->end(),res);
2701 DataArrayDouble *MEDCouplingFieldDiscretizationKriging::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfTargetPoints) const
2703 if(!arr || !arr->isAllocated())
2704 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationKriging::getValueOnMulti : input array is null or not allocated !");
2705 int nbOfRows(getNumberOfMeshPlaces(mesh));
2706 if(arr->getNumberOfTuples()!=nbOfRows)
2708 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationKriging::getValueOnMulti : input array does not have correct number of tuples ! Excepted " << nbOfRows << " having " << arr->getNumberOfTuples() << " !";
2709 throw INTERP_KERNEL::Exception(oss.str().c_str());
2711 int nbCols(-1),nbCompo(arr->getNumberOfComponents());
2712 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> m(computeEvaluationMatrixOnGivenPts(mesh,loc,nbOfTargetPoints,nbCols));
2713 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret(DataArrayDouble::New());
2714 ret->alloc(nbOfTargetPoints,nbCompo);
2715 INTERP_KERNEL::matrixProduct(m->begin(),nbOfTargetPoints,nbCols,arr->begin(),nbOfRows,nbCompo,ret->getPointer());
2719 void MEDCouplingFieldDiscretizationKriging::reprQuickOverview(std::ostream& stream) const
2721 stream << "Kriging spatial discretization.";
2725 * Returns the matrix of size nbRows = \a nbOfTargetPoints and \a nbCols = \a nbCols. This matrix is useful if
2727 * \return the new result matrix to be deallocated by the caller.
2729 DataArrayDouble *MEDCouplingFieldDiscretizationKriging::computeEvaluationMatrixOnGivenPts(const MEDCouplingMesh *mesh, const double *loc, int nbOfTargetPoints, int& nbCols) const
2731 int isDrift(-1),nbRows(-1);
2732 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrixInv(computeInverseMatrix(mesh,isDrift,nbRows));
2734 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coords=getLocalizationOfDiscValues(mesh);
2735 int nbOfPts(coords->getNumberOfTuples()),dimension(coords->getNumberOfComponents());
2736 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> locArr=DataArrayDouble::New();
2737 locArr->useArray(loc,false,CPP_DEALLOC,nbOfTargetPoints,dimension);
2740 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrix2=coords->buildEuclidianDistanceDenseMatrixWith(locArr);
2741 operateOnDenseMatrix(mesh->getSpaceDimension(),nbOfTargetPoints*nbOfPts,matrix2->getPointer());
2743 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrix3=DataArrayDouble::New();
2744 matrix3->alloc(nbOfTargetPoints*nbRows,1);
2745 double *work=matrix3->getPointer();
2746 const double *workCst(matrix2->begin()),*workCst2(loc);
2747 for(int i=0;i<nbOfTargetPoints;i++,workCst+=nbOfPts,workCst2+=isDrift-1)
2749 for(int j=0;j<nbOfPts;j++)
2750 work[i*nbRows+j]=workCst[j];
2751 work[i*nbRows+nbOfPts]=1.0;
2752 for(int j=0;j<isDrift-1;j++)
2753 work[i*nbRows+(nbOfPts+1+j)]=workCst2[j];
2755 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret(DataArrayDouble::New());
2756 ret->alloc(nbOfTargetPoints,nbRows);
2757 INTERP_KERNEL::matrixProduct(matrix3->begin(),nbOfTargetPoints,nbRows,matrixInv->begin(),nbRows,nbRows,ret->getPointer());
2758 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret2(DataArrayDouble::New());
2759 ret2->alloc(nbOfTargetPoints*nbOfPts,1);
2760 workCst=ret->begin(); work=ret2->getPointer();
2761 for(int i=0;i<nbOfTargetPoints;i++,workCst+=nbRows)
2762 work=std::copy(workCst,workCst+nbOfPts,work);
2767 * This method returns the square matrix of size \a matSz that is the inverse of the kriging matrix. The returned matrix can returned all the coeffs of kriging
2768 * when multiplied by the vector of values attached to each point.
2770 * \param [out] isDrift return if drift coefficients are present in the returned vector of coefficients. If different from 0 there is presence of drift coefficients.
2771 * \param [out] matSz the size of returned square matrix
2772 * \return the new result matrix to be deallocated by the caller.
2774 DataArrayDouble *MEDCouplingFieldDiscretizationKriging::computeInverseMatrix(const MEDCouplingMesh *mesh, int& isDrift, int& matSz) const
2777 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationKriging::computeVectorOfCoefficients : NULL input mesh !");
2778 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coords=getLocalizationOfDiscValues(mesh);
2779 int nbOfPts=coords->getNumberOfTuples();
2780 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrix=coords->buildEuclidianDistanceDenseMatrix();
2781 operateOnDenseMatrix(mesh->getSpaceDimension(),nbOfPts*nbOfPts,matrix->getPointer());
2783 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrixWithDrift=performDrift(matrix,coords,isDrift);
2784 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrixInv=DataArrayDouble::New();
2785 matSz=nbOfPts+isDrift;
2786 matrixInv->alloc(matSz*matSz,1);
2787 INTERP_KERNEL::inverseMatrix(matrixWithDrift->getConstPointer(),matSz,matrixInv->getPointer());
2788 return matrixInv.retn();
2792 * 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
2793 * number of tuples should be equal to the number of representing points in \a mesh.
2795 * \param [in] mesh is the sources of nodes on which kriging will be done regarding the parameters and the value of \c this->getSpaceDimension()
2796 * \param [in] arr input field DataArrayDouble whose number of tuples must be equal to the number of nodes in \a mesh
2797 * \param [out] isDrift return if drift coefficients are present in the returned vector of coefficients. If different from 0 there is presence of drift coefficients.
2798 * Whatever the value of \a isDrift the number of tuples of returned DataArrayDouble will be equal to \c arr->getNumberOfTuples() + \a isDrift.
2799 * \return a newly allocated array containing coefficients including or not drift coefficient at the end depending the value of \a isDrift parameter.
2801 DataArrayDouble *MEDCouplingFieldDiscretizationKriging::computeVectorOfCoefficients(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, int& isDrift) const
2804 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrixInv(computeInverseMatrix(mesh,isDrift,nbRows));
2805 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> KnewiK=DataArrayDouble::New();
2806 KnewiK->alloc(nbRows*1,1);
2807 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arr2=DataArrayDouble::New();
2808 arr2->alloc(nbRows*1,1);
2809 double *work=std::copy(arr->begin(),arr->end(),arr2->getPointer());
2810 std::fill(work,work+isDrift,0.);
2811 INTERP_KERNEL::matrixProduct(matrixInv->getConstPointer(),nbRows,nbRows,arr2->getConstPointer(),nbRows,1,KnewiK->getPointer());
2812 return KnewiK.retn();
2816 * 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.
2818 * \param [in] spaceDimension space dimension of the input mesh on which the Kriging has to be performed
2819 * \param [in] nbOfElems is the result of the product of nb of rows and the nb of columns of matrix \a matrixPtr
2820 * \param [in,out] matrixPtr is the dense matrix whose on each values the operation will be applied
2822 void MEDCouplingFieldDiscretizationKriging::operateOnDenseMatrix(int spaceDimension, int nbOfElems, double *matrixPtr) const
2824 switch(spaceDimension)
2828 for(int i=0;i<nbOfElems;i++)
2830 double val=matrixPtr[i];
2831 matrixPtr[i]=val*val*val;
2837 for(int i=0;i<nbOfElems;i++)
2839 double val=matrixPtr[i];
2841 matrixPtr[i]=val*val*log(val);
2847 //nothing here : it is not a bug g(h)=h with spaceDim 3.
2851 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationKriging::operateOnDenseMatrix : only dimension 1, 2 and 3 implemented !");
2856 * Starting from a square matrix \a matr, this method returns a newly allocated dense square matrix whose \a matr is included in returned matrix
2857 * in the top left corner, and in the remaining returned matrix the parameters to take into account about the kriging drift.
2858 * For the moment only linear srift is implemented.
2860 * \param [in] arr the position of points were input mesh geometry is considered for Kriging
2861 * \param [in] matr input matrix whose drift part will be added
2862 * \param [out] delta the difference between the size of the output matrix and the input matrix \a matr.
2863 * \return a newly allocated matrix bigger than input matrix \a matr.
2865 DataArrayDouble *MEDCouplingFieldDiscretizationKriging::performDrift(const DataArrayDouble *matr, const DataArrayDouble *arr, int& delta) const
2867 int spaceDimension=arr->getNumberOfComponents();
2868 delta=spaceDimension+1;
2869 int szOfMatrix=arr->getNumberOfTuples();
2870 if(szOfMatrix*szOfMatrix!=matr->getNumberOfTuples())
2871 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationKriging::performDrift : invalid size");
2872 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
2873 ret->alloc((szOfMatrix+delta)*(szOfMatrix+delta),1);
2874 const double *srcWork=matr->getConstPointer();
2875 const double *srcWork2=arr->getConstPointer();
2876 double *destWork=ret->getPointer();
2877 for(int i=0;i<szOfMatrix;i++)
2879 destWork=std::copy(srcWork,srcWork+szOfMatrix,destWork);
2880 srcWork+=szOfMatrix;
2882 destWork=std::copy(srcWork2,srcWork2+spaceDimension,destWork);
2883 srcWork2+=spaceDimension;
2885 std::fill(destWork,destWork+szOfMatrix,1.); destWork+=szOfMatrix;
2886 std::fill(destWork,destWork+spaceDimension+1,0.); destWork+=spaceDimension+1;
2887 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arrNoI=arr->toNoInterlace();
2888 srcWork2=arrNoI->getConstPointer();
2889 for(int i=0;i<spaceDimension;i++)
2891 destWork=std::copy(srcWork2,srcWork2+szOfMatrix,destWork);
2892 srcWork2+=szOfMatrix;
2893 std::fill(destWork,destWork+spaceDimension+1,0.);
2894 destWork+=spaceDimension+1;