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_POINT1[1]={0.};
69 const double MEDCouplingFieldDiscretizationGaussNE::FGP_SEG2[2]={1.,1.};
70 const double MEDCouplingFieldDiscretizationGaussNE::FGP_SEG3[3]={0.5555555555555556,0.8888888888888888,0.5555555555555556};
71 const double MEDCouplingFieldDiscretizationGaussNE::FGP_SEG4[4]={0.347854845137454,0.347854845137454,0.652145154862546,0.652145154862546};
72 const double MEDCouplingFieldDiscretizationGaussNE::FGP_TRI3[3]={0.16666666666666666,0.16666666666666666,0.16666666666666666};
73 const double MEDCouplingFieldDiscretizationGaussNE::FGP_TRI6[6]={0.0549758718227661,0.0549758718227661,0.0549758718227661,0.11169079483905,0.11169079483905,0.11169079483905};
74 const double MEDCouplingFieldDiscretizationGaussNE::FGP_TRI7[7]={0.062969590272413,0.062969590272413,0.062969590272413,0.066197076394253,0.066197076394253,0.066197076394253,0.1125};
75 const double MEDCouplingFieldDiscretizationGaussNE::FGP_QUAD4[4]={1.,1.,1.,1.};
76 const double MEDCouplingFieldDiscretizationGaussNE::FGP_QUAD8[8]={1.,1.,1.,1.,1.,1.,1.,1.};
77 const double MEDCouplingFieldDiscretizationGaussNE::FGP_QUAD9[9]={0.30864197530864196,0.30864197530864196,0.30864197530864196,0.30864197530864196,0.49382716049382713,0.49382716049382713,0.49382716049382713,0.49382716049382713,0.7901234567901234};
78 const double MEDCouplingFieldDiscretizationGaussNE::FGP_TETRA4[4]={0.041666666666666664,0.041666666666666664,0.041666666666666664,0.041666666666666664};
79 const double MEDCouplingFieldDiscretizationGaussNE::FGP_TETRA10[10]={1.,1.,1.,1.,1.,1.,1.,1.,1.,1.};//to check
80 const double MEDCouplingFieldDiscretizationGaussNE::FGP_PENTA6[6]={0.16666666666666666,0.16666666666666666,0.16666666666666666,0.16666666666666666,0.16666666666666666,0.16666666666666666};
81 const double MEDCouplingFieldDiscretizationGaussNE::FGP_PENTA15[15]={1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.};//to check
82 const double MEDCouplingFieldDiscretizationGaussNE::FGP_HEXA8[8]={1.,1.,1.,1.,1.,1.,1.,1.};
83 const double MEDCouplingFieldDiscretizationGaussNE::FGP_HEXA20[20]={1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.};
84 const double MEDCouplingFieldDiscretizationGaussNE::FGP_HEXA27[27]={1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.};
85 const double MEDCouplingFieldDiscretizationGaussNE::FGP_PYRA5[5]={0.13333333333333333,0.13333333333333333,0.13333333333333333,0.13333333333333333,0.13333333333333333};
86 const double MEDCouplingFieldDiscretizationGaussNE::FGP_PYRA13[13]={1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.};//to check
87 const double MEDCouplingFieldDiscretizationGaussNE::REF_SEG2[2]={-1.,1.};
88 const double MEDCouplingFieldDiscretizationGaussNE::REF_SEG3[3]={-1.,1.,0.};
89 const double MEDCouplingFieldDiscretizationGaussNE::REF_SEG4[4]={-1.,1.,-0.3333333333333333,0.3333333333333333};
90 const double MEDCouplingFieldDiscretizationGaussNE::REF_TRI3[6]={0.,0.,1.,0.,0.,1.};
91 const double MEDCouplingFieldDiscretizationGaussNE::REF_TRI6[12]={0.,0.,1.,0.,0.,1.,0.5,0.,0.5,0.5,0.,0.5};
92 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};
93 const double MEDCouplingFieldDiscretizationGaussNE::REF_QUAD4[8]={-1.,-1.,1.,-1.,1.,1.,-1.,1.};
94 const double MEDCouplingFieldDiscretizationGaussNE::REF_QUAD8[16]={-1.,-1.,1.,-1.,1.,1.,-1.,1.,0.,-1.,1.,0.,0.,1.,-1.,0.};
95 const double MEDCouplingFieldDiscretizationGaussNE::REF_QUAD9[18]={-1.,-1.,1.,-1.,1.,1.,-1.,1.,0.,-1.,1.,0.,0.,1.,-1.,0.,0.,0.};
96 const double MEDCouplingFieldDiscretizationGaussNE::REF_TETRA4[12]={0.,1.,0.,0.,0.,1.,0.,0.,0.,1.,0.,0.};
97 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.};
98 const double MEDCouplingFieldDiscretizationGaussNE::REF_PENTA6[18]={-1.,1.,0.,-1.,0.,1.,-1.,0.,0.,1.,1.,0.,1.,0.,1.,1.,0.,0.};
99 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.};
100 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.};
101 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.};
102 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.,-1.,0.,-1.,0.,1.,-1.,1.,0.,-1.,0.,-1.,-1.,-1.,0.,1.,0.,1.,1.,1.,0.,1.,0.,-1.,1.,-1.,-1.,0.,-1.,1.,0.,1.,1.,0.,1.,-1.,0.,0.,0.,-1.,-1.,0.,0.,0.,1.,0.,1.,0.,0.,0.,-1.,0.,0.,0.,1.,0.,0.,0.};
103 const double MEDCouplingFieldDiscretizationGaussNE::REF_PYRA5[15]={1.,0.,0.,0.,1.,0.,-1.,0.,0.,0.,-1.,0.,0.,0.,1.};
104 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};
105 const double MEDCouplingFieldDiscretizationGaussNE::LOC_SEG2[2]={0.577350269189626,-0.577350269189626};
106 const double MEDCouplingFieldDiscretizationGaussNE::LOC_SEG3[3]={-0.774596669241,0.,0.774596669241};
107 const double MEDCouplingFieldDiscretizationGaussNE::LOC_SEG4[4]={0.339981043584856,-0.339981043584856,0.861136311594053,-0.861136311594053};
108 const double MEDCouplingFieldDiscretizationGaussNE::LOC_TRI3[6]={0.16666666666666667,0.16666666666666667,0.6666666666666667,0.16666666666666667,0.16666666666666667,0.6666666666666667};
109 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};
110 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};
111 const double MEDCouplingFieldDiscretizationGaussNE::LOC_QUAD4[8]={-0.774596669241483,-0.774596669241483,0.774596669241483,-0.774596669241483,0.774596669241483,0.774596669241483,-0.774596669241483,0.774596669241483};
112 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.};
113 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.};
114 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};
115 const double MEDCouplingFieldDiscretizationGaussNE::LOC_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.};//to check
116 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.};
117 const double MEDCouplingFieldDiscretizationGaussNE::LOC_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.};//to check
118 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};
119 const double MEDCouplingFieldDiscretizationGaussNE::LOC_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.};//to check
120 const double MEDCouplingFieldDiscretizationGaussNE::LOC_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.,-1.,0.,-1.,0.,1.,-1.,1.,0.,-1.,0.,-1.,-1.,-1.,0.,1.,0.,1.,1.,1.,0.,1.,0.,-1.,1.,-1.,-1.,0.,-1.,1.,0.,1.,1.,0.,1.,-1.,0.,0.,0.,-1.,-1.,0.,0.,0.,1.,0.,1.,0.,0.,0.,-1.,0.,0.,0.,1.,0.,0.,0.};
121 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};
122 const double MEDCouplingFieldDiscretizationGaussNE::LOC_PYRA13[39]={1.,0.,0.,0.,-1.,0.,-1.,0.,0.,0.,1.,0.,0.,0.,0.999999999999,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};//to check 0.99999... to avoid nan ! on node #4 of PYRA13
124 MEDCouplingFieldDiscretization::MEDCouplingFieldDiscretization():_precision(DFLT_PRECISION)
128 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretization::New(TypeOfField type)
132 case MEDCouplingFieldDiscretizationP0::TYPE:
133 return new MEDCouplingFieldDiscretizationP0;
134 case MEDCouplingFieldDiscretizationP1::TYPE:
135 return new MEDCouplingFieldDiscretizationP1;
136 case MEDCouplingFieldDiscretizationGauss::TYPE:
137 return new MEDCouplingFieldDiscretizationGauss;
138 case MEDCouplingFieldDiscretizationGaussNE::TYPE:
139 return new MEDCouplingFieldDiscretizationGaussNE;
140 case MEDCouplingFieldDiscretizationKriging::TYPE:
141 return new MEDCouplingFieldDiscretizationKriging;
143 throw INTERP_KERNEL::Exception("Choosen discretization is not implemented yet.");
147 TypeOfField MEDCouplingFieldDiscretization::GetTypeOfFieldFromStringRepr(const std::string& repr)
149 if(repr==MEDCouplingFieldDiscretizationP0::REPR)
150 return MEDCouplingFieldDiscretizationP0::TYPE;
151 if(repr==MEDCouplingFieldDiscretizationP1::REPR)
152 return MEDCouplingFieldDiscretizationP1::TYPE;
153 if(repr==MEDCouplingFieldDiscretizationGauss::REPR)
154 return MEDCouplingFieldDiscretizationGauss::TYPE;
155 if(repr==MEDCouplingFieldDiscretizationGaussNE::REPR)
156 return MEDCouplingFieldDiscretizationGaussNE::TYPE;
157 if(repr==MEDCouplingFieldDiscretizationKriging::REPR)
158 return MEDCouplingFieldDiscretizationKriging::TYPE;
159 throw INTERP_KERNEL::Exception("Representation does not match with any field discretization !");
162 bool MEDCouplingFieldDiscretization::isEqual(const MEDCouplingFieldDiscretization *other, double eps) const
165 return isEqualIfNotWhy(other,eps,reason);
168 bool MEDCouplingFieldDiscretization::isEqualWithoutConsideringStr(const MEDCouplingFieldDiscretization *other, double eps) const
170 return isEqual(other,eps);
174 * This method is an alias of MEDCouplingFieldDiscretization::clone. It is only here for coherency with all the remaining of MEDCoupling.
175 * \sa MEDCouplingFieldDiscretization::clone.
177 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretization::deepCpy() const
183 * For all field discretization excepted GaussPts the [ \a startCellIds, \a endCellIds ) has no impact on the cloned instance.
185 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretization::clonePart(const int *startCellIds, const int *endCellIds) const
191 * For all field discretization excepted GaussPts the slice( \a beginCellId, \a endCellIds, \a stepCellId ) has no impact on the cloned instance.
193 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretization::clonePartRange(int beginCellIds, int endCellIds, int stepCellIds) const
199 * Excepted for MEDCouplingFieldDiscretizationPerCell no underlying TimeLabel object : nothing to do in generally.
201 void MEDCouplingFieldDiscretization::updateTime() const
205 std::size_t MEDCouplingFieldDiscretization::getHeapMemorySizeWithoutChildren() const
210 std::vector<const BigMemoryObject *> MEDCouplingFieldDiscretization::getDirectChildren() const
212 return std::vector<const BigMemoryObject *>();
216 * Computes normL1 of DataArrayDouble instance arr.
217 * @param res output parameter expected to be of size arr->getNumberOfComponents();
218 * @throw when the field discretization fails on getMeasure fields (gauss points for example)
220 void MEDCouplingFieldDiscretization::normL1(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, double *res) const
222 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=getMeasureField(mesh,true);
223 int nbOfCompo=arr->getNumberOfComponents();
224 int nbOfElems=getNumberOfTuples(mesh);
225 std::fill(res,res+nbOfCompo,0.);
226 const double *arrPtr=arr->getConstPointer();
227 const double *volPtr=vol->getArray()->getConstPointer();
229 for(int i=0;i<nbOfElems;i++)
231 double v=fabs(volPtr[i]);
232 for(int j=0;j<nbOfCompo;j++)
233 res[j]+=fabs(arrPtr[i*nbOfCompo+j])*v;
236 std::transform(res,res+nbOfCompo,res,std::bind2nd(std::multiplies<double>(),1./deno));
240 * Computes normL2 of DataArrayDouble instance arr.
241 * @param res output parameter expected to be of size arr->getNumberOfComponents();
242 * @throw when the field discretization fails on getMeasure fields (gauss points for example)
244 void MEDCouplingFieldDiscretization::normL2(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, double *res) const
246 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=getMeasureField(mesh,true);
247 int nbOfCompo=arr->getNumberOfComponents();
248 int nbOfElems=getNumberOfTuples(mesh);
249 std::fill(res,res+nbOfCompo,0.);
250 const double *arrPtr=arr->getConstPointer();
251 const double *volPtr=vol->getArray()->getConstPointer();
253 for(int i=0;i<nbOfElems;i++)
255 double v=fabs(volPtr[i]);
256 for(int j=0;j<nbOfCompo;j++)
257 res[j]+=arrPtr[i*nbOfCompo+j]*arrPtr[i*nbOfCompo+j]*v;
260 std::transform(res,res+nbOfCompo,res,std::bind2nd(std::multiplies<double>(),1./deno));
261 std::transform(res,res+nbOfCompo,res,std::ptr_fun<double,double>(std::sqrt));
265 * Computes integral of DataArrayDouble instance arr.
266 * @param res output parameter expected to be of size arr->getNumberOfComponents();
267 * @throw when the field discretization fails on getMeasure fields (gauss points for example)
269 void MEDCouplingFieldDiscretization::integral(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, bool isWAbs, double *res) const
272 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretization::integral : mesh is NULL !");
274 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretization::integral : input array is NULL !");
275 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=getMeasureField(mesh,isWAbs);
276 int nbOfCompo=arr->getNumberOfComponents();
277 int nbOfElems=getNumberOfTuples(mesh);
278 if(nbOfElems!=arr->getNumberOfTuples())
280 std::ostringstream oss; oss << "MEDCouplingFieldDiscretization::integral : field is not correct ! number of tuples in array is " << arr->getNumberOfTuples();
281 oss << " whereas number of tuples expected is " << nbOfElems << " !";
282 throw INTERP_KERNEL::Exception(oss.str().c_str());
284 std::fill(res,res+nbOfCompo,0.);
285 const double *arrPtr=arr->getConstPointer();
286 const double *volPtr=vol->getArray()->getConstPointer();
287 INTERP_KERNEL::AutoPtr<double> tmp=new double[nbOfCompo];
288 for (int i=0;i<nbOfElems;i++)
290 std::transform(arrPtr+i*nbOfCompo,arrPtr+(i+1)*nbOfCompo,(double *)tmp,std::bind2nd(std::multiplies<double>(),volPtr[i]));
291 std::transform((double *)tmp,(double *)tmp+nbOfCompo,res,res,std::plus<double>());
296 * This method is strictly equivalent to MEDCouplingFieldDiscretization::buildSubMeshData except that it is optimized for input defined as a range of cell ids.
298 * \param [out] beginOut Valid only if \a di is NULL
299 * \param [out] endOut Valid only if \a di is NULL
300 * \param [out] stepOut Valid only if \a di is NULL
301 * \param [out] di is an array returned that specifies entity ids (nodes, cells, Gauss points... ) in array if no output range is foundable.
303 * \sa MEDCouplingFieldDiscretization::buildSubMeshData
305 MEDCouplingMesh *MEDCouplingFieldDiscretization::buildSubMeshDataRange(const MEDCouplingMesh *mesh, int beginCellIds, int endCellIds, int stepCellIds, int& beginOut, int& endOut, int& stepOut, DataArrayInt *&di) const
307 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> da=DataArrayInt::Range(beginCellIds,endCellIds,stepCellIds);
308 return buildSubMeshData(mesh,da->begin(),da->end(),di);
311 void MEDCouplingFieldDiscretization::getSerializationIntArray(DataArrayInt *& arr) const
319 void MEDCouplingFieldDiscretization::getTinySerializationIntInformation(std::vector<int>& tinyInfo) const
326 void MEDCouplingFieldDiscretization::getTinySerializationDbleInformation(std::vector<double>& tinyInfo) const
330 void MEDCouplingFieldDiscretization::resizeForUnserialization(const std::vector<int>& tinyInfo, DataArrayInt *& arr)
338 void MEDCouplingFieldDiscretization::finishUnserialization(const std::vector<double>& tinyInfo)
343 * This method is typically the first step of renumbering. The implementation is empty it is not a bug only gauss is impacted
344 * virtualy by this method.
346 void MEDCouplingFieldDiscretization::renumberCells(const int *old2NewBg, bool check)
350 double MEDCouplingFieldDiscretization::getIJK(const MEDCouplingMesh *mesh, const DataArrayDouble *da, int cellId, int nodeIdInCell, int compoId) const
352 throw INTERP_KERNEL::Exception("getIJK Invalid ! only for GaussPoint and GaussNE discretizations !");
355 void MEDCouplingFieldDiscretization::setGaussLocalizationOnType(const MEDCouplingMesh *m, INTERP_KERNEL::NormalizedCellType type, const std::vector<double>& refCoo,
356 const std::vector<double>& gsCoo, const std::vector<double>& wg)
358 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
361 void MEDCouplingFieldDiscretization::setGaussLocalizationOnCells(const MEDCouplingMesh *m, const int *begin, const int *end, const std::vector<double>& refCoo,
362 const std::vector<double>& gsCoo, const std::vector<double>& wg)
364 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
367 void MEDCouplingFieldDiscretization::clearGaussLocalizations()
369 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
372 MEDCouplingGaussLocalization& MEDCouplingFieldDiscretization::getGaussLocalization(int locId)
374 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
377 const MEDCouplingGaussLocalization& MEDCouplingFieldDiscretization::getGaussLocalization(int locId) const
379 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
382 int MEDCouplingFieldDiscretization::getNbOfGaussLocalization() const
384 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
387 int MEDCouplingFieldDiscretization::getGaussLocalizationIdOfOneCell(int cellId) const
389 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
392 int MEDCouplingFieldDiscretization::getGaussLocalizationIdOfOneType(INTERP_KERNEL::NormalizedCellType type) const
394 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
397 std::set<int> MEDCouplingFieldDiscretization::getGaussLocalizationIdsOfOneType(INTERP_KERNEL::NormalizedCellType type) const
399 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
402 void MEDCouplingFieldDiscretization::getCellIdsHavingGaussLocalization(int locId, std::vector<int>& cellIds) const
404 throw INTERP_KERNEL::Exception("Invalid method for the corresponding field discretization : available only for GaussPoint discretization !");
407 void MEDCouplingFieldDiscretization::RenumberEntitiesFromO2NArr(double eps, const int *old2NewPtr, int newNbOfEntity, DataArrayDouble *arr, const std::string& msg)
410 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretization::RenumberEntitiesFromO2NArr : input array is NULL !");
411 int oldNbOfElems=arr->getNumberOfTuples();
412 int nbOfComp=arr->getNumberOfComponents();
413 int newNbOfTuples=newNbOfEntity;
414 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arrCpy=arr->deepCpy();
415 const double *ptSrc=arrCpy->getConstPointer();
416 arr->reAlloc(newNbOfTuples);
417 double *ptToFill=arr->getPointer();
418 std::fill(ptToFill,ptToFill+nbOfComp*newNbOfTuples,std::numeric_limits<double>::max());
419 INTERP_KERNEL::AutoPtr<double> tmp=new double[nbOfComp];
420 for(int i=0;i<oldNbOfElems;i++)
422 int newNb=old2NewPtr[i];
423 if(newNb>=0)//if newNb<0 the node is considered as out.
425 if(std::find_if(ptToFill+newNb*nbOfComp,ptToFill+(newNb+1)*nbOfComp,std::bind2nd(std::not_equal_to<double>(),std::numeric_limits<double>::max()))
426 ==ptToFill+(newNb+1)*nbOfComp)
427 std::copy(ptSrc+i*nbOfComp,ptSrc+(i+1)*nbOfComp,ptToFill+newNb*nbOfComp);
430 std::transform(ptSrc+i*nbOfComp,ptSrc+(i+1)*nbOfComp,ptToFill+newNb*nbOfComp,(double *)tmp,std::minus<double>());
431 std::transform((double *)tmp,((double *)tmp)+nbOfComp,(double *)tmp,std::ptr_fun<double,double>(fabs));
432 //if(!std::equal(ptSrc+i*nbOfComp,ptSrc+(i+1)*nbOfComp,ptToFill+newNb*nbOfComp))
433 if(*std::max_element((double *)tmp,((double *)tmp)+nbOfComp)>eps)
435 std::ostringstream oss;
436 oss << msg << " " << i << " and " << std::find(old2NewPtr,old2NewPtr+i,newNb)-old2NewPtr
437 << " have been merged and " << msg << " field on them are different !";
438 throw INTERP_KERNEL::Exception(oss.str().c_str());
445 void MEDCouplingFieldDiscretization::RenumberEntitiesFromN2OArr(const int *new2OldPtr, int new2OldSz, DataArrayDouble *arr, const std::string& msg)
447 int nbOfComp=arr->getNumberOfComponents();
448 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arrCpy=arr->deepCpy();
449 const double *ptSrc=arrCpy->getConstPointer();
450 arr->reAlloc(new2OldSz);
451 double *ptToFill=arr->getPointer();
452 for(int i=0;i<new2OldSz;i++)
454 int oldNb=new2OldPtr[i];
455 std::copy(ptSrc+oldNb*nbOfComp,ptSrc+(oldNb+1)*nbOfComp,ptToFill+i*nbOfComp);
459 MEDCouplingFieldDiscretization::~MEDCouplingFieldDiscretization()
463 TypeOfField MEDCouplingFieldDiscretizationP0::getEnum() const
469 * This method is simply called by MEDCouplingFieldDiscretization::deepCpy. It performs the deep copy of \a this.
471 * \sa MEDCouplingFieldDiscretization::deepCpy.
473 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationP0::clone() const
475 return new MEDCouplingFieldDiscretizationP0;
478 std::string MEDCouplingFieldDiscretizationP0::getStringRepr() const
480 return std::string(REPR);
483 const char *MEDCouplingFieldDiscretizationP0::getRepr() const
488 bool MEDCouplingFieldDiscretizationP0::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
492 reason="other spatial discretization is NULL, and this spatial discretization (P0) is defined.";
495 const MEDCouplingFieldDiscretizationP0 *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationP0 *>(other);
498 reason="Spatial discrtization of this is ON_CELLS, which is not the case of other.";
502 int MEDCouplingFieldDiscretizationP0::getNumberOfTuples(const MEDCouplingMesh *mesh) const
505 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getNumberOfTuples : NULL input mesh !");
506 return mesh->getNumberOfCells();
510 * This method returns the number of tuples regarding exclusively the input code \b without \b using \b a \b mesh \b in \b input.
511 * The input code coherency is also checked regarding spatial discretization of \a this.
512 * If an incoherency is detected, an exception will be thrown. If the input code is coherent, the number of tuples expected is returned.
513 * 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).
515 int MEDCouplingFieldDiscretizationP0::getNumberOfTuplesExpectedRegardingCode(const std::vector<int>& code, const std::vector<const DataArrayInt *>& idsPerType) const
518 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getNumberOfTuplesExpectedRegardingCode : invalid input code !");
519 int nbOfSplit=(int)idsPerType.size();
520 int nbOfTypes=(int)code.size()/3;
522 for(int i=0;i<nbOfTypes;i++)
524 int nbOfEltInChunk=code[3*i+1];
526 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getNumberOfTuplesExpectedRegardingCode : invalid input code ! presence of negative value in a type !");
530 if(pos<0 || pos>=nbOfSplit)
532 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationP0::getNumberOfTuplesExpectedRegardingCode : input code points to pos " << pos << " in typeid " << i << " ! Should be in [0," << nbOfSplit << ") !";
533 throw INTERP_KERNEL::Exception(oss.str().c_str());
535 const DataArrayInt *ids(idsPerType[pos]);
536 if(!ids || !ids->isAllocated() || ids->getNumberOfComponents()!=1 || ids->getNumberOfTuples()!=nbOfEltInChunk || ids->getMinValueInArray()<0)
538 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationP0::getNumberOfTuplesExpectedRegardingCode : input pfl chunck at pos " << pos << " should have " << i << " tuples and one component and with ids all >=0 !";
539 throw INTERP_KERNEL::Exception(oss.str().c_str());
547 int MEDCouplingFieldDiscretizationP0::getNumberOfMeshPlaces(const MEDCouplingMesh *mesh) const
550 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getNumberOfMeshPlaces : NULL input mesh !");
551 return mesh->getNumberOfCells();
554 DataArrayInt *MEDCouplingFieldDiscretizationP0::getOffsetArr(const MEDCouplingMesh *mesh) const
557 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getOffsetArr : NULL input mesh !");
558 int nbOfTuples=mesh->getNumberOfCells();
559 DataArrayInt *ret=DataArrayInt::New();
560 ret->alloc(nbOfTuples+1,1);
565 void MEDCouplingFieldDiscretizationP0::renumberArraysForCell(const MEDCouplingMesh *mesh, const std::vector<DataArray *>& arrays,
566 const int *old2NewBg, bool check)
569 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::renumberArraysForCell : NULL input mesh !");
570 const int *array=old2NewBg;
572 array=DataArrayInt::CheckAndPreparePermutation(old2NewBg,old2NewBg+mesh->getNumberOfCells());
573 for(std::vector<DataArray *>::const_iterator it=arrays.begin();it!=arrays.end();it++)
576 (*it)->renumberInPlace(array);
579 free(const_cast<int *>(array));
582 DataArrayDouble *MEDCouplingFieldDiscretizationP0::getLocalizationOfDiscValues(const MEDCouplingMesh *mesh) const
585 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getLocalizationOfDiscValues : NULL input mesh !");
586 return mesh->getBarycenterAndOwner();
589 void MEDCouplingFieldDiscretizationP0::computeMeshRestrictionFromTupleIds(const MEDCouplingMesh *mesh, const int *tupleIdsBg, const int *tupleIdsEnd,
590 DataArrayInt *&cellRestriction, DataArrayInt *&trueTupleRestriction) const
593 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::computeMeshRestrictionFromTupleIds : NULL input mesh !");
594 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> tmp=DataArrayInt::New();
595 tmp->alloc((int)std::distance(tupleIdsBg,tupleIdsEnd),1);
596 std::copy(tupleIdsBg,tupleIdsEnd,tmp->getPointer());
597 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> tmp2(tmp->deepCpy());
598 cellRestriction=tmp.retn();
599 trueTupleRestriction=tmp2.retn();
602 void MEDCouplingFieldDiscretizationP0::reprQuickOverview(std::ostream& stream) const
604 stream << "P0 spatial discretization.";
607 void MEDCouplingFieldDiscretizationP0::checkCompatibilityWithNature(NatureOfField nat) const
611 void MEDCouplingFieldDiscretizationP0::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArray *da) const
614 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::checkCoherencyBetween : NULL input mesh or DataArray !");
615 if(mesh->getNumberOfCells()!=da->getNumberOfTuples())
617 std::ostringstream message;
618 message << "Field on cells invalid because there are " << mesh->getNumberOfCells();
619 message << " cells in mesh and " << da->getNumberOfTuples() << " tuples in field !";
620 throw INTERP_KERNEL::Exception(message.str().c_str());
624 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationP0::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
627 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getMeasureField : mesh instance specified is NULL !");
628 return mesh->getMeasureField(isAbs);
631 void MEDCouplingFieldDiscretizationP0::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
634 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getValueOn : NULL input mesh !");
635 int id=mesh->getCellContainingPoint(loc,_precision);
637 throw INTERP_KERNEL::Exception("Specified point is detected outside of mesh : unable to apply P0::getValueOn !");
638 arr->getTuple(id,res);
641 void MEDCouplingFieldDiscretizationP0::getValueOnPos(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, int i, int j, int k, double *res) const
643 const MEDCouplingCMesh *meshC=dynamic_cast<const MEDCouplingCMesh *>(mesh);
645 throw INTERP_KERNEL::Exception("P0::getValueOnPos is only accessible for structured meshes !");
646 int id=meshC->getCellIdFromPos(i,j,k);
647 arr->getTuple(id,res);
650 DataArrayDouble *MEDCouplingFieldDiscretizationP0::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const
653 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getValueOnMulti : NULL input mesh !");
654 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> eltsArr,eltsIndexArr;
655 mesh->getCellsContainingPoints(loc,nbOfPoints,_precision,eltsArr,eltsIndexArr);
656 const int *elts(eltsArr->begin()),*eltsIndex(eltsIndexArr->begin());
657 int spaceDim=mesh->getSpaceDimension();
658 int nbOfComponents=arr->getNumberOfComponents();
659 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
660 ret->alloc(nbOfPoints,nbOfComponents);
661 double *ptToFill=ret->getPointer();
662 for(int i=0;i<nbOfPoints;i++,ptToFill+=nbOfComponents)
663 if(eltsIndex[i+1]-eltsIndex[i]>=1)
664 arr->getTuple(elts[eltsIndex[i]],ptToFill);
667 std::ostringstream oss; oss << "Point #" << i << " with coordinates : (";
668 std::copy(loc+i*spaceDim,loc+(i+1)*spaceDim,std::ostream_iterator<double>(oss,", "));
669 oss << ") detected outside mesh : unable to apply P0::getValueOnMulti ! ";
670 throw INTERP_KERNEL::Exception(oss.str().c_str());
676 * Nothing to do. It's not a bug.
678 void MEDCouplingFieldDiscretizationP0::renumberValuesOnNodes(double , const int *, int newNbOfNodes, DataArrayDouble *) const
682 void MEDCouplingFieldDiscretizationP0::renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const
684 RenumberEntitiesFromO2NArr(epsOnVals,old2New,newSz,arr,"Cell");
687 void MEDCouplingFieldDiscretizationP0::renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const
689 RenumberEntitiesFromN2OArr(new2old,newSz,arr,"Cell");
693 * This method returns a tuple ids selection from cell ids selection [start;end).
694 * This method is called by MEDCouplingFieldDiscretizationP0::buildSubMeshData to return parameter \b di.
695 * Here for P0 it's very simple !
697 * \return a newly allocated array containing ids to select into the DataArrayDouble of the field.
700 DataArrayInt *MEDCouplingFieldDiscretizationP0::computeTupleIdsToSelectFromCellIds(const MEDCouplingMesh *mesh, const int *startCellIds, const int *endCellIds) const
702 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret=DataArrayInt::New();
703 ret->alloc((int)std::distance(startCellIds,endCellIds),1);
704 std::copy(startCellIds,endCellIds,ret->getPointer());
709 * This method returns a submesh of 'mesh' instance constituting cell ids contained in array defined as an interval [start;end).
710 * @param di is an array returned that specifies entity ids (here cells ids) in mesh 'mesh' of entity in returned submesh.
711 * Example : The first cell id of returned mesh has the (*di)[0] id in 'mesh'
713 * \sa MEDCouplingFieldDiscretizationP0::buildSubMeshDataRange
715 MEDCouplingMesh *MEDCouplingFieldDiscretizationP0::buildSubMeshData(const MEDCouplingMesh *mesh, const int *start, const int *end, DataArrayInt *&di) const
718 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::buildSubMeshData : NULL input mesh !");
719 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPart(start,end);
720 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> diSafe=DataArrayInt::New();
721 diSafe->alloc((int)std::distance(start,end),1);
722 std::copy(start,end,diSafe->getPointer());
728 * This method is strictly equivalent to MEDCouplingFieldDiscretizationP0::buildSubMeshData except that it is optimized for input defined as a range of cell ids.
730 * \param [out] beginOut Valid only if \a di is NULL
731 * \param [out] endOut Valid only if \a di is NULL
732 * \param [out] stepOut Valid only if \a di is NULL
733 * \param [out] di is an array returned that specifies entity ids (nodes, cells, Gauss points... ) in array if no output range is foundable.
735 * \sa MEDCouplingFieldDiscretizationP0::buildSubMeshData
737 MEDCouplingMesh *MEDCouplingFieldDiscretizationP0::buildSubMeshDataRange(const MEDCouplingMesh *mesh, int beginCellIds, int endCellIds, int stepCellIds, int& beginOut, int& endOut, int& stepOut, DataArrayInt *&di) const
740 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::buildSubMeshDataRange : NULL input mesh !");
741 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPartRange(beginCellIds,endCellIds,stepCellIds);
742 di=0; beginOut=beginCellIds; endOut=endCellIds; stepOut=stepCellIds;
746 int MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuples(const MEDCouplingMesh *mesh) const
749 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::getNumberOfTuples : NULL input mesh !");
750 return mesh->getNumberOfNodes();
754 * This method returns the number of tuples regarding exclusively the input code \b without \b using \b a \b mesh \b in \b input.
755 * The input code coherency is also checked regarding spatial discretization of \a this.
756 * If an incoherency is detected, an exception will be thrown. If the input code is coherent, the number of tuples expected is returned.
757 * 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).
759 int MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuplesExpectedRegardingCode(const std::vector<int>& code, const std::vector<const DataArrayInt *>& idsPerType) const
762 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuplesExpectedRegardingCode : invalid input code !");
763 int nbOfSplit=(int)idsPerType.size();
764 int nbOfTypes=(int)code.size()/3;
766 for(int i=0;i<nbOfTypes;i++)
768 int nbOfEltInChunk=code[3*i+1];
770 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuplesExpectedRegardingCode : invalid input code ! presence of negative value in a type !");
774 if(pos<0 || pos>=nbOfSplit)
776 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuplesExpectedRegardingCode : input code points to pos " << pos << " in typeid " << i << " ! Should be in [0," << nbOfSplit << ") !";
777 throw INTERP_KERNEL::Exception(oss.str().c_str());
779 const DataArrayInt *ids(idsPerType[pos]);
780 if(!ids || !ids->isAllocated() || ids->getNumberOfComponents()!=1 || ids->getNumberOfTuples()!=nbOfEltInChunk || ids->getMinValueInArray()<0)
782 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuplesExpectedRegardingCode : input pfl chunck at pos " << pos << " should have " << i << " tuples and one component and with ids all >=0 !";
783 throw INTERP_KERNEL::Exception(oss.str().c_str());
791 int MEDCouplingFieldDiscretizationOnNodes::getNumberOfMeshPlaces(const MEDCouplingMesh *mesh) const
794 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::getNumberOfMeshPlaces : NULL input mesh !");
795 return mesh->getNumberOfNodes();
799 * Nothing to do here.
801 void MEDCouplingFieldDiscretizationOnNodes::renumberArraysForCell(const MEDCouplingMesh *, const std::vector<DataArray *>& arrays,
802 const int *old2NewBg, bool check)
806 DataArrayInt *MEDCouplingFieldDiscretizationOnNodes::getOffsetArr(const MEDCouplingMesh *mesh) const
809 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::getOffsetArr : NULL input mesh !");
810 int nbOfTuples=mesh->getNumberOfNodes();
811 DataArrayInt *ret=DataArrayInt::New();
812 ret->alloc(nbOfTuples+1,1);
817 DataArrayDouble *MEDCouplingFieldDiscretizationOnNodes::getLocalizationOfDiscValues(const MEDCouplingMesh *mesh) const
820 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::getLocalizationOfDiscValues : NULL input mesh !");
821 return mesh->getCoordinatesAndOwner();
824 void MEDCouplingFieldDiscretizationOnNodes::computeMeshRestrictionFromTupleIds(const MEDCouplingMesh *mesh, const int *tupleIdsBg, const int *tupleIdsEnd,
825 DataArrayInt *&cellRestriction, DataArrayInt *&trueTupleRestriction) const
828 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationOnNodes::computeMeshRestrictionFromTupleIds : NULL input mesh !");
829 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret1=mesh->getCellIdsFullyIncludedInNodeIds(tupleIdsBg,tupleIdsEnd);
830 const MEDCouplingUMesh *meshc=dynamic_cast<const MEDCouplingUMesh *>(mesh);
832 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationOnNodes::computeMeshRestrictionFromTupleIds : trying to subpart field on nodes by node ids ! Your mesh has to be unstructured !");
833 MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> meshPart=static_cast<MEDCouplingUMesh *>(meshc->buildPartOfMySelf(ret1->begin(),ret1->end(),true));
834 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret2=meshPart->computeFetchedNodeIds();
835 cellRestriction=ret1.retn();
836 trueTupleRestriction=ret2.retn();
839 void MEDCouplingFieldDiscretizationOnNodes::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArray *da) const
842 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::checkCoherencyBetween : NULL input mesh or DataArray !");
843 if(mesh->getNumberOfNodes()!=da->getNumberOfTuples())
845 std::ostringstream message;
846 message << "Field on nodes invalid because there are " << mesh->getNumberOfNodes();
847 message << " nodes in mesh and " << da->getNumberOfTuples() << " tuples in field !";
848 throw INTERP_KERNEL::Exception(message.str().c_str());
853 * This method returns a submesh of 'mesh' instance constituting cell ids contained in array defined as an interval [start;end).
854 * @param di is an array returned that specifies entity ids (here nodes ids) in mesh 'mesh' of entity in returned submesh.
855 * Example : The first node id of returned mesh has the (*di)[0] id in 'mesh'
857 MEDCouplingMesh *MEDCouplingFieldDiscretizationOnNodes::buildSubMeshData(const MEDCouplingMesh *mesh, const int *start, const int *end, DataArrayInt *&di) const
860 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationNodes::buildSubMeshData : NULL input mesh !");
861 DataArrayInt *diTmp=0;
862 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPartAndReduceNodes(start,end,diTmp);
863 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> diTmpSafe(diTmp);
864 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> di2=diTmpSafe->invertArrayO2N2N2O(ret->getNumberOfNodes());
870 * This method is strictly equivalent to MEDCouplingFieldDiscretizationNodes::buildSubMeshData except that it is optimized for input defined as a range of cell ids.
872 * \param [out] beginOut Valid only if \a di is NULL
873 * \param [out] endOut Valid only if \a di is NULL
874 * \param [out] stepOut Valid only if \a di is NULL
875 * \param [out] di is an array returned that specifies entity ids (nodes, cells, Gauss points... ) in array if no output range is foundable.
877 * \sa MEDCouplingFieldDiscretizationNodes::buildSubMeshData
879 MEDCouplingMesh *MEDCouplingFieldDiscretizationOnNodes::buildSubMeshDataRange(const MEDCouplingMesh *mesh, int beginCellIds, int endCellIds, int stepCellIds, int& beginOut, int& endOut, int& stepOut, DataArrayInt *&di) const
882 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationOnNodes::buildSubMeshDataRange : NULL input mesh !");
883 DataArrayInt *diTmp=0;
884 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPartRangeAndReduceNodes(beginCellIds,endCellIds,stepCellIds,beginOut,endOut,stepOut,diTmp);
887 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> diTmpSafe(diTmp);
888 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> di2=diTmpSafe->invertArrayO2N2N2O(ret->getNumberOfNodes());
895 * This method returns a tuple ids selection from cell ids selection [start;end).
896 * This method is called by MEDCouplingFieldDiscretizationOnNodes::buildSubMeshData to return parameter \b di.
897 * Here for P1 only nodes fetched by submesh of mesh[startCellIds:endCellIds) is returned !
899 * \return a newly allocated array containing ids to select into the DataArrayDouble of the field.
902 DataArrayInt *MEDCouplingFieldDiscretizationOnNodes::computeTupleIdsToSelectFromCellIds(const MEDCouplingMesh *mesh, const int *startCellIds, const int *endCellIds) const
905 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP1::computeTupleIdsToSelectFromCellIds : NULL input mesh !");
906 const MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> umesh=mesh->buildUnstructured();
907 MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> umesh2=static_cast<MEDCouplingUMesh *>(umesh->buildPartOfMySelf(startCellIds,endCellIds,true));
908 return umesh2->computeFetchedNodeIds();
911 void MEDCouplingFieldDiscretizationOnNodes::renumberValuesOnNodes(double epsOnVals, const int *old2NewPtr, int newNbOfNodes, DataArrayDouble *arr) const
913 RenumberEntitiesFromO2NArr(epsOnVals,old2NewPtr,newNbOfNodes,arr,"Node");
917 * Nothing to do it's not a bug.
919 void MEDCouplingFieldDiscretizationOnNodes::renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const
924 * Nothing to do it's not a bug.
926 void MEDCouplingFieldDiscretizationOnNodes::renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const
930 void MEDCouplingFieldDiscretizationOnNodes::getValueOnPos(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, int i, int j, int k, double *res) const
932 const MEDCouplingCMesh *meshC=dynamic_cast<const MEDCouplingCMesh *>(mesh);
934 throw INTERP_KERNEL::Exception("OnNodes::getValueOnPos(i,j,k) is only accessible for structured meshes !");
935 int id=meshC->getNodeIdFromPos(i,j,k);
936 arr->getTuple(id,res);
939 TypeOfField MEDCouplingFieldDiscretizationP1::getEnum() const
945 * This method is simply called by MEDCouplingFieldDiscretization::deepCpy. It performs the deep copy of \a this.
947 * \sa MEDCouplingFieldDiscretization::deepCpy.
949 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationP1::clone() const
951 return new MEDCouplingFieldDiscretizationP1;
954 std::string MEDCouplingFieldDiscretizationP1::getStringRepr() const
956 return std::string(REPR);
959 const char *MEDCouplingFieldDiscretizationP1::getRepr() const
964 bool MEDCouplingFieldDiscretizationP1::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
968 reason="other spatial discretization is NULL, and this spatial discretization (P1) is defined.";
971 const MEDCouplingFieldDiscretizationP1 *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationP1 *>(other);
974 reason="Spatial discrtization of this is ON_NODES, which is not the case of other.";
978 void MEDCouplingFieldDiscretizationP1::checkCompatibilityWithNature(NatureOfField nat) const
980 if(nat!=ConservativeVolumic)
981 throw INTERP_KERNEL::Exception("Invalid nature for P1 field : expected ConservativeVolumic !");
984 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationP1::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
987 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP1::getMeasureField : mesh instance specified is NULL !");
988 return mesh->getMeasureFieldOnNode(isAbs);
991 void MEDCouplingFieldDiscretizationP1::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
994 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP1::getValueOn : NULL input mesh !");
995 int id=mesh->getCellContainingPoint(loc,_precision);
997 throw INTERP_KERNEL::Exception("Specified point is detected outside of mesh : unable to apply P1::getValueOn !");
998 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(id);
999 if(type!=INTERP_KERNEL::NORM_SEG2 && type!=INTERP_KERNEL::NORM_TRI3 && type!=INTERP_KERNEL::NORM_TETRA4)
1000 throw INTERP_KERNEL::Exception("P1 getValueOn is not specified for not simplex cells !");
1001 getValueInCell(mesh,id,arr,loc,res);
1005 * This method localizes a point defined by 'loc' in a cell with id 'cellId' into mesh 'mesh'.
1006 * The result is put into res expected to be of size at least arr->getNumberOfComponents()
1008 void MEDCouplingFieldDiscretizationP1::getValueInCell(const MEDCouplingMesh *mesh, int cellId, const DataArrayDouble *arr, const double *loc, double *res) const
1011 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP1::getValueInCell : NULL input mesh !");
1012 std::vector<int> conn;
1013 std::vector<double> coo;
1014 mesh->getNodeIdsOfCell(cellId,conn);
1015 for(std::vector<int>::const_iterator iter=conn.begin();iter!=conn.end();iter++)
1016 mesh->getCoordinatesOfNode(*iter,coo);
1017 int spaceDim=mesh->getSpaceDimension();
1018 std::size_t nbOfNodes=conn.size();
1019 std::vector<const double *> vec(nbOfNodes);
1020 for(std::size_t i=0;i<nbOfNodes;i++)
1021 vec[i]=&coo[i*spaceDim];
1022 INTERP_KERNEL::AutoPtr<double> tmp=new double[nbOfNodes];
1023 INTERP_KERNEL::barycentric_coords(vec,loc,tmp);
1024 int sz=arr->getNumberOfComponents();
1025 INTERP_KERNEL::AutoPtr<double> tmp2=new double[sz];
1026 std::fill(res,res+sz,0.);
1027 for(std::size_t i=0;i<nbOfNodes;i++)
1029 arr->getTuple(conn[i],(double *)tmp2);
1030 std::transform((double *)tmp2,((double *)tmp2)+sz,(double *)tmp2,std::bind2nd(std::multiplies<double>(),tmp[i]));
1031 std::transform(res,res+sz,(double *)tmp2,res,std::plus<double>());
1035 DataArrayDouble *MEDCouplingFieldDiscretizationP1::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const
1038 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP1::getValueOnMulti : NULL input mesh !");
1039 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> eltsArr,eltsIndexArr;
1040 mesh->getCellsContainingPoints(loc,nbOfPoints,_precision,eltsArr,eltsIndexArr);
1041 const int *elts(eltsArr->begin()),*eltsIndex(eltsIndexArr->begin());
1042 int spaceDim=mesh->getSpaceDimension();
1043 int nbOfComponents=arr->getNumberOfComponents();
1044 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
1045 ret->alloc(nbOfPoints,nbOfComponents);
1046 double *ptToFill=ret->getPointer();
1047 for(int i=0;i<nbOfPoints;i++)
1048 if(eltsIndex[i+1]-eltsIndex[i]>=1)
1049 getValueInCell(mesh,elts[eltsIndex[i]],arr,loc+i*spaceDim,ptToFill+i*nbOfComponents);
1052 std::ostringstream oss; oss << "Point #" << i << " with coordinates : (";
1053 std::copy(loc+i*spaceDim,loc+(i+1)*spaceDim,std::ostream_iterator<double>(oss,", "));
1054 oss << ") detected outside mesh : unable to apply P1::getValueOnMulti ! ";
1055 throw INTERP_KERNEL::Exception(oss.str().c_str());
1060 void MEDCouplingFieldDiscretizationP1::reprQuickOverview(std::ostream& stream) const
1062 stream << "P1 spatial discretization.";
1065 MEDCouplingFieldDiscretizationPerCell::MEDCouplingFieldDiscretizationPerCell():_discr_per_cell(0)
1069 MEDCouplingFieldDiscretizationPerCell::~MEDCouplingFieldDiscretizationPerCell()
1072 _discr_per_cell->decrRef();
1076 * This constructor deep copies ParaMEDMEM::DataArrayInt instance from other (if any).
1078 MEDCouplingFieldDiscretizationPerCell::MEDCouplingFieldDiscretizationPerCell(const MEDCouplingFieldDiscretizationPerCell& other, const int *startCellIds, const int *endCellIds):_discr_per_cell(0)
1080 DataArrayInt *arr=other._discr_per_cell;
1083 if(startCellIds==0 && endCellIds==0)
1084 _discr_per_cell=arr->deepCpy();
1086 _discr_per_cell=arr->selectByTupleIdSafe(startCellIds,endCellIds);
1090 MEDCouplingFieldDiscretizationPerCell::MEDCouplingFieldDiscretizationPerCell(const MEDCouplingFieldDiscretizationPerCell& other, int beginCellIds, int endCellIds, int stepCellIds):_discr_per_cell(0)
1092 DataArrayInt *arr=other._discr_per_cell;
1095 _discr_per_cell=arr->selectByTupleId2(beginCellIds,endCellIds,stepCellIds);
1099 void MEDCouplingFieldDiscretizationPerCell::updateTime() const
1102 updateTimeWith(*_discr_per_cell);
1105 std::size_t MEDCouplingFieldDiscretizationPerCell::getHeapMemorySizeWithoutChildren() const
1107 std::size_t ret(MEDCouplingFieldDiscretization::getHeapMemorySizeWithoutChildren());
1111 std::vector<const BigMemoryObject *> MEDCouplingFieldDiscretizationPerCell::getDirectChildren() const
1113 std::vector<const BigMemoryObject *> ret(MEDCouplingFieldDiscretization::getDirectChildren());
1115 ret.push_back(_discr_per_cell);
1119 void MEDCouplingFieldDiscretizationPerCell::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArray *da) const
1121 if(!_discr_per_cell)
1122 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell has no discretization per cell !");
1124 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell::checkCoherencyBetween : NULL input mesh or DataArray !");
1125 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1126 if(nbOfTuples!=mesh->getNumberOfCells())
1127 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell has a discretization per cell but it's not matching the underlying mesh !");
1130 bool MEDCouplingFieldDiscretizationPerCell::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
1134 reason="other spatial discretization is NULL, and this spatial discretization (PerCell) is defined.";
1137 const MEDCouplingFieldDiscretizationPerCell *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationPerCell *>(other);
1140 reason="Spatial discretization of this is ON_GAUSS, which is not the case of other.";
1143 if(_discr_per_cell==0)
1144 return otherC->_discr_per_cell==0;
1145 if(otherC->_discr_per_cell==0)
1147 bool ret=_discr_per_cell->isEqualIfNotWhy(*otherC->_discr_per_cell,reason);
1149 reason.insert(0,"Field discretization per cell DataArrayInt given the discid per cell :");
1153 bool MEDCouplingFieldDiscretizationPerCell::isEqualWithoutConsideringStr(const MEDCouplingFieldDiscretization *other, double eps) const
1155 const MEDCouplingFieldDiscretizationPerCell *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationPerCell *>(other);
1158 if(_discr_per_cell==0)
1159 return otherC->_discr_per_cell==0;
1160 if(otherC->_discr_per_cell==0)
1162 return _discr_per_cell->isEqualWithoutConsideringStr(*otherC->_discr_per_cell);
1166 * This method is typically the first step of renumbering. The impact on _discr_per_cell is necessary here.
1167 * virtualy by this method.
1169 void MEDCouplingFieldDiscretizationPerCell::renumberCells(const int *old2NewBg, bool check)
1171 int nbCells=_discr_per_cell->getNumberOfTuples();
1172 const int *array=old2NewBg;
1174 array=DataArrayInt::CheckAndPreparePermutation(old2NewBg,old2NewBg+nbCells);
1176 DataArrayInt *dpc=_discr_per_cell->renumber(array);
1177 _discr_per_cell->decrRef();
1178 _discr_per_cell=dpc;
1181 free(const_cast<int *>(array));
1184 void MEDCouplingFieldDiscretizationPerCell::buildDiscrPerCellIfNecessary(const MEDCouplingMesh *mesh)
1187 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell::buildDiscrPerCellIfNecessary : NULL input mesh !");
1188 if(!_discr_per_cell)
1190 _discr_per_cell=DataArrayInt::New();
1191 int nbTuples=mesh->getNumberOfCells();
1192 _discr_per_cell->alloc(nbTuples,1);
1193 int *ptr=_discr_per_cell->getPointer();
1194 std::fill(ptr,ptr+nbTuples,DFT_INVALID_LOCID_VALUE);
1198 void MEDCouplingFieldDiscretizationPerCell::checkNoOrphanCells() const
1200 if(!_discr_per_cell)
1201 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell::checkNoOrphanCells : no discretization defined !");
1202 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> test=_discr_per_cell->getIdsEqual(DFT_INVALID_LOCID_VALUE);
1203 if(test->getNumberOfTuples()!=0)
1204 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell::checkNoOrphanCells : presence of orphan cells !");
1208 * This method is useful when 'this' describes a field discretization with several gauss discretization on a \b same cell type.
1209 * For example same NORM_TRI3 cells having 6 gauss points and others with 12 gauss points.
1210 * This method returns 2 arrays with same size : the return value and 'locIds' output parameter.
1211 * 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.
1212 * The return vector contains a set of newly created instance to deal with.
1213 * The returned vector represents a \b partition of cells ids with a gauss discretization set.
1215 * If no descretization is set in 'this' and exception will be thrown.
1217 std::vector<DataArrayInt *> MEDCouplingFieldDiscretizationPerCell::splitIntoSingleGaussDicrPerCellType(std::vector<int>& locIds) const
1219 if(!_discr_per_cell)
1220 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationPerCell::splitIntoSingleGaussDicrPerCellType : no descretization set !");
1221 return _discr_per_cell->partitionByDifferentValues(locIds);
1224 const DataArrayInt *MEDCouplingFieldDiscretizationPerCell::getArrayOfDiscIds() const
1226 return _discr_per_cell;
1229 void MEDCouplingFieldDiscretizationPerCell::setArrayOfDiscIds(const DataArrayInt *adids)
1231 if(adids!=_discr_per_cell)
1234 _discr_per_cell->decrRef();
1235 _discr_per_cell=const_cast<DataArrayInt *>(adids);
1237 _discr_per_cell->incrRef();
1242 MEDCouplingFieldDiscretizationGauss::MEDCouplingFieldDiscretizationGauss()
1246 MEDCouplingFieldDiscretizationGauss::MEDCouplingFieldDiscretizationGauss(const MEDCouplingFieldDiscretizationGauss& other, const int *startCellIds, const int *endCellIds):MEDCouplingFieldDiscretizationPerCell(other,startCellIds,endCellIds),_loc(other._loc)
1250 MEDCouplingFieldDiscretizationGauss::MEDCouplingFieldDiscretizationGauss(const MEDCouplingFieldDiscretizationGauss& other, int beginCellIds, int endCellIds, int stepCellIds):MEDCouplingFieldDiscretizationPerCell(other,beginCellIds,endCellIds,stepCellIds),_loc(other._loc)
1254 TypeOfField MEDCouplingFieldDiscretizationGauss::getEnum() const
1259 bool MEDCouplingFieldDiscretizationGauss::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
1263 reason="other spatial discretization is NULL, and this spatial discretization (Gauss) is defined.";
1266 const MEDCouplingFieldDiscretizationGauss *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationGauss *>(other);
1269 reason="Spatial discrtization of this is ON_GAUSS, which is not the case of other.";
1272 if(!MEDCouplingFieldDiscretizationPerCell::isEqualIfNotWhy(other,eps,reason))
1274 if(_loc.size()!=otherC->_loc.size())
1276 reason="Gauss spatial discretization : localization sizes differ";
1279 std::size_t sz=_loc.size();
1280 for(std::size_t i=0;i<sz;i++)
1281 if(!_loc[i].isEqual(otherC->_loc[i],eps))
1283 std::ostringstream oss; oss << "Gauss spatial discretization : Localization #" << i << " differ from this to other.";
1290 bool MEDCouplingFieldDiscretizationGauss::isEqualWithoutConsideringStr(const MEDCouplingFieldDiscretization *other, double eps) const
1292 const MEDCouplingFieldDiscretizationGauss *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationGauss *>(other);
1295 if(!MEDCouplingFieldDiscretizationPerCell::isEqualWithoutConsideringStr(other,eps))
1297 if(_loc.size()!=otherC->_loc.size())
1299 std::size_t sz=_loc.size();
1300 for(std::size_t i=0;i<sz;i++)
1301 if(!_loc[i].isEqual(otherC->_loc[i],eps))
1307 * This method is simply called by MEDCouplingFieldDiscretization::deepCpy. It performs the deep copy of \a this.
1309 * \sa MEDCouplingFieldDiscretization::deepCpy.
1311 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationGauss::clone() const
1313 return new MEDCouplingFieldDiscretizationGauss(*this);
1316 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationGauss::clonePart(const int *startCellIds, const int *endCellIds) const
1318 return new MEDCouplingFieldDiscretizationGauss(*this,startCellIds,endCellIds);
1321 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationGauss::clonePartRange(int beginCellIds, int endCellIds, int stepCellIds) const
1323 return new MEDCouplingFieldDiscretizationGauss(*this,beginCellIds,endCellIds,stepCellIds);
1326 std::string MEDCouplingFieldDiscretizationGauss::getStringRepr() const
1328 std::ostringstream oss; oss << REPR << "." << std::endl;
1331 if(_discr_per_cell->isAllocated())
1333 oss << "Discretization per cell : ";
1334 std::copy(_discr_per_cell->begin(),_discr_per_cell->end(),std::ostream_iterator<int>(oss,", "));
1338 oss << "Presence of " << _loc.size() << " localizations." << std::endl;
1340 for(std::vector<MEDCouplingGaussLocalization>::const_iterator it=_loc.begin();it!=_loc.end();it++,i++)
1342 oss << "+++++ Localization #" << i << " +++++" << std::endl;
1343 oss << (*it).getStringRepr();
1344 oss << "++++++++++" << std::endl;
1349 std::size_t MEDCouplingFieldDiscretizationGauss::getHeapMemorySizeWithoutChildren() const
1351 std::size_t ret(MEDCouplingFieldDiscretizationPerCell::getHeapMemorySizeWithoutChildren());
1352 ret+=_loc.capacity()*sizeof(MEDCouplingGaussLocalization);
1353 for(std::vector<MEDCouplingGaussLocalization>::const_iterator it=_loc.begin();it!=_loc.end();it++)
1354 ret+=(*it).getMemorySize();
1358 const char *MEDCouplingFieldDiscretizationGauss::getRepr() const
1364 * This method returns the number of tuples regarding exclusively the input code \b without \b using \b a \b mesh \b in \b input.
1365 * The input code coherency is also checked regarding spatial discretization of \a this.
1366 * If an incoherency is detected, an exception will be thrown. If the input code is coherent, the number of tuples expected is returned.
1367 * 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).
1369 int MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode(const std::vector<int>& code, const std::vector<const DataArrayInt *>& idsPerType) const
1371 if(!_discr_per_cell || !_discr_per_cell->isAllocated() || _discr_per_cell->getNumberOfComponents()!=1)
1372 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode");
1373 if(code.size()%3!=0)
1374 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode : invalid input code !");
1375 int nbOfSplit=(int)idsPerType.size();
1376 int nbOfTypes=(int)code.size()/3;
1378 for(int i=0;i<nbOfTypes;i++)
1380 int nbOfEltInChunk=code[3*i+1];
1381 if(nbOfEltInChunk<0)
1382 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode : invalid input code ! presence of negative value in a type !");
1383 int pos=code[3*i+2];
1386 if(pos<0 || pos>=nbOfSplit)
1388 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode : input code points to pos " << pos << " in typeid " << i << " ! Should be in [0," << nbOfSplit << ") !";
1389 throw INTERP_KERNEL::Exception(oss.str().c_str());
1391 const DataArrayInt *ids(idsPerType[pos]);
1392 if(!ids || !ids->isAllocated() || ids->getNumberOfComponents()!=1 || ids->getNumberOfTuples()!=nbOfEltInChunk || ids->getMinValueInArray()<0)
1394 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode : input pfl chunck at pos " << pos << " should have " << i << " tuples and one component and with ids all >=0 !";
1395 throw INTERP_KERNEL::Exception(oss.str().c_str());
1398 ret+=nbOfEltInChunk;
1400 if(ret!=_discr_per_cell->getNumberOfTuples())
1402 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getNumberOfTuplesExpectedRegardingCode : input code points to " << ret << " cells whereas discretization percell array lgth is " << _discr_per_cell->getNumberOfTuples() << " !";
1403 throw INTERP_KERNEL::Exception(oss.str().c_str());
1405 return getNumberOfTuples(0);//0 is not an error ! It is to be sure that input mesh is not used
1408 int MEDCouplingFieldDiscretizationGauss::getNumberOfTuples(const MEDCouplingMesh *) const
1411 if (_discr_per_cell == 0)
1412 throw INTERP_KERNEL::Exception("Discretization is not initialized!");
1413 const int *dcPtr=_discr_per_cell->getConstPointer();
1414 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1415 int maxSz=(int)_loc.size();
1416 for(const int *w=dcPtr;w!=dcPtr+nbOfTuples;w++)
1418 if(*w>=0 && *w<maxSz)
1419 ret+=_loc[*w].getNumberOfGaussPt();
1422 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getNumberOfTuples : At cell #" << std::distance(dcPtr,w) << " localization id is " << *w << " should be in [0," << maxSz << ") !";
1423 throw INTERP_KERNEL::Exception(oss.str().c_str());
1429 int MEDCouplingFieldDiscretizationGauss::getNumberOfMeshPlaces(const MEDCouplingMesh *mesh) const
1432 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getNumberOfMeshPlaces : NULL input mesh !");
1433 return mesh->getNumberOfCells();
1437 * This method is redevelopped for performance reasons, but it is equivalent to a call to MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField
1438 * and a call to DataArrayDouble::computeOffsets2 on the returned array.
1440 DataArrayInt *MEDCouplingFieldDiscretizationGauss::getOffsetArr(const MEDCouplingMesh *mesh) const
1443 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getOffsetArr : NULL input mesh !");
1444 int nbOfTuples=mesh->getNumberOfCells();
1445 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret=DataArrayInt::New();
1446 ret->alloc(nbOfTuples+1,1);
1447 int *retPtr=ret->getPointer();
1448 const int *start=_discr_per_cell->getConstPointer();
1449 if(_discr_per_cell->getNumberOfTuples()!=nbOfTuples)
1450 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getOffsetArr : mismatch between the mesh and the discretization ids array length !");
1451 int maxPossible=(int)_loc.size();
1453 for(int i=0;i<nbOfTuples;i++,start++)
1455 if(*start>=0 && *start<maxPossible)
1456 retPtr[i+1]=retPtr[i]+_loc[*start].getNumberOfGaussPt();
1459 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getOffsetArr : At position #" << i << " the locid = " << *start << " whereas it should be in [0," << maxPossible << ") !";
1460 throw INTERP_KERNEL::Exception(oss.str().c_str());
1466 void MEDCouplingFieldDiscretizationGauss::renumberArraysForCell(const MEDCouplingMesh *mesh, const std::vector<DataArray *>& arrays,
1467 const int *old2NewBg, bool check)
1470 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::renumberArraysForCell : NULL input mesh !");
1471 const int *array=old2NewBg;
1473 array=DataArrayInt::CheckAndPreparePermutation(old2NewBg,old2NewBg+mesh->getNumberOfCells());
1474 int nbOfCells=_discr_per_cell->getNumberOfTuples();
1475 int nbOfTuples=getNumberOfTuples(0);
1476 const int *dcPtr=_discr_per_cell->getConstPointer();
1477 int *array2=new int[nbOfTuples];//stores the final conversion array old2New to give to arrays in renumberInPlace.
1478 int *array3=new int[nbOfCells];//store for each cell in present dcp array (already renumbered) the offset needed by each cell in new numbering.
1480 for(int i=1;i<nbOfCells;i++)
1481 array3[i]=array3[i-1]+_loc[dcPtr[i-1]].getNumberOfGaussPt();
1483 for(int i=0;i<nbOfCells;i++)
1485 int nbOfGaussPt=_loc[dcPtr[array[i]]].getNumberOfGaussPt();
1486 for(int k=0;k<nbOfGaussPt;k++,j++)
1487 array2[j]=array3[array[i]]+k;
1490 for(std::vector<DataArray *>::const_iterator it=arrays.begin();it!=arrays.end();it++)
1492 (*it)->renumberInPlace(array2);
1495 free(const_cast<int*>(array));
1498 DataArrayDouble *MEDCouplingFieldDiscretizationGauss::getLocalizationOfDiscValues(const MEDCouplingMesh *mesh) const
1501 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getLocalizationOfDiscValues : NULL input mesh !");
1502 checkNoOrphanCells();
1503 MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> umesh=mesh->buildUnstructured();//in general do nothing
1504 int nbOfTuples=getNumberOfTuples(mesh);
1505 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
1506 int spaceDim=mesh->getSpaceDimension();
1507 ret->alloc(nbOfTuples,spaceDim);
1508 std::vector< int > locIds;
1509 std::vector<DataArrayInt *> parts=splitIntoSingleGaussDicrPerCellType(locIds);
1510 std::vector< MEDCouplingAutoRefCountObjectPtr<DataArrayInt> > parts2(parts.size());
1511 std::copy(parts.begin(),parts.end(),parts2.begin());
1512 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> offsets=buildNbOfGaussPointPerCellField();
1513 offsets->computeOffsets();
1514 const int *ptrOffsets=offsets->getConstPointer();
1515 const double *coords=umesh->getCoords()->getConstPointer();
1516 const int *connI=umesh->getNodalConnectivityIndex()->getConstPointer();
1517 const int *conn=umesh->getNodalConnectivity()->getConstPointer();
1518 double *valsToFill=ret->getPointer();
1519 for(std::size_t i=0;i<parts2.size();i++)
1521 INTERP_KERNEL::GaussCoords calculator;
1523 const MEDCouplingGaussLocalization& cli=_loc[locIds[i]];//curLocInfo
1524 INTERP_KERNEL::NormalizedCellType typ=cli.getType();
1525 const std::vector<double>& wg=cli.getWeights();
1526 calculator.addGaussInfo(typ,INTERP_KERNEL::CellModel::GetCellModel(typ).getDimension(),
1527 &cli.getGaussCoords()[0],(int)wg.size(),&cli.getRefCoords()[0],
1528 INTERP_KERNEL::CellModel::GetCellModel(typ).getNumberOfNodes());
1530 int nbt=parts2[i]->getNumberOfTuples();
1531 for(const int *w=parts2[i]->getConstPointer();w!=parts2[i]->getConstPointer()+nbt;w++)
1532 calculator.calculateCoords(cli.getType(),coords,spaceDim,conn+connI[*w]+1,valsToFill+spaceDim*(ptrOffsets[*w]));
1534 ret->copyStringInfoFrom(*umesh->getCoords());
1538 void MEDCouplingFieldDiscretizationGauss::computeMeshRestrictionFromTupleIds(const MEDCouplingMesh *mesh, const int *tupleIdsBg, const int *tupleIdsEnd,
1539 DataArrayInt *&cellRestriction, DataArrayInt *&trueTupleRestriction) const
1542 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::computeMeshRestrictionFromTupleIds : NULL input mesh !");
1543 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> tmp=DataArrayInt::New(); tmp->alloc((int)std::distance(tupleIdsBg,tupleIdsEnd),1);
1544 std::copy(tupleIdsBg,tupleIdsEnd,tmp->getPointer());
1546 tmp=tmp->buildUnique();
1547 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=buildNbOfGaussPointPerCellField();
1548 nbOfNodesPerCell->computeOffsets2();
1549 nbOfNodesPerCell->searchRangesInListOfIds(tmp,cellRestriction,trueTupleRestriction);
1555 void MEDCouplingFieldDiscretizationGauss::checkCompatibilityWithNature(NatureOfField nat) const
1559 void MEDCouplingFieldDiscretizationGauss::getTinySerializationIntInformation(std::vector<int>& tinyInfo) const
1563 val=_discr_per_cell->getNumberOfTuples();
1564 tinyInfo.push_back(val);
1565 tinyInfo.push_back((int)_loc.size());
1567 tinyInfo.push_back(-1);
1569 tinyInfo.push_back(_loc[0].getDimension());
1570 for(std::vector<MEDCouplingGaussLocalization>::const_iterator iter=_loc.begin();iter!=_loc.end();iter++)
1571 (*iter).pushTinySerializationIntInfo(tinyInfo);
1574 void MEDCouplingFieldDiscretizationGauss::getTinySerializationDbleInformation(std::vector<double>& tinyInfo) const
1576 for(std::vector<MEDCouplingGaussLocalization>::const_iterator iter=_loc.begin();iter!=_loc.end();iter++)
1577 (*iter).pushTinySerializationDblInfo(tinyInfo);
1580 void MEDCouplingFieldDiscretizationGauss::getSerializationIntArray(DataArrayInt *& arr) const
1584 arr=_discr_per_cell;
1587 void MEDCouplingFieldDiscretizationGauss::resizeForUnserialization(const std::vector<int>& tinyInfo, DataArrayInt *& arr)
1589 int val=tinyInfo[0];
1592 _discr_per_cell=DataArrayInt::New();
1593 _discr_per_cell->alloc(val,1);
1597 arr=_discr_per_cell;
1598 int nbOfLoc=tinyInfo[1];
1600 int dim=tinyInfo[2];
1603 delta=((int)tinyInfo.size()-3)/nbOfLoc;
1604 for(int i=0;i<nbOfLoc;i++)
1606 std::vector<int> tmp(tinyInfo.begin()+3+i*delta,tinyInfo.begin()+3+(i+1)*delta);
1607 MEDCouplingGaussLocalization elt=MEDCouplingGaussLocalization::BuildNewInstanceFromTinyInfo(dim,tmp);
1608 _loc.push_back(elt);
1612 void MEDCouplingFieldDiscretizationGauss::finishUnserialization(const std::vector<double>& tinyInfo)
1614 double *tmp=new double[tinyInfo.size()];
1615 std::copy(tinyInfo.begin(),tinyInfo.end(),tmp);
1616 const double *work=tmp;
1617 for(std::vector<MEDCouplingGaussLocalization>::iterator iter=_loc.begin();iter!=_loc.end();iter++)
1618 work=(*iter).fillWithValues(work);
1622 double MEDCouplingFieldDiscretizationGauss::getIJK(const MEDCouplingMesh *mesh, const DataArrayDouble *da, int cellId, int nodeIdInCell, int compoId) const
1624 int offset=getOffsetOfCell(cellId);
1625 return da->getIJ(offset+nodeIdInCell,compoId);
1628 void MEDCouplingFieldDiscretizationGauss::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArray *da) const
1631 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::checkCoherencyBetween : NULL input mesh or DataArray !");
1632 MEDCouplingFieldDiscretizationPerCell::checkCoherencyBetween(mesh,da);
1633 for(std::vector<MEDCouplingGaussLocalization>::const_iterator iter=_loc.begin();iter!=_loc.end();iter++)
1634 (*iter).checkCoherency();
1635 int nbOfDesc=(int)_loc.size();
1636 int nbOfCells=mesh->getNumberOfCells();
1637 const int *dc=_discr_per_cell->getConstPointer();
1638 for(int i=0;i<nbOfCells;i++)
1642 std::ostringstream oss; oss << "Cell # " << i << " of mesh \"" << mesh->getName() << "\" has an undefined gauss location ! Should never happend !";
1643 throw INTERP_KERNEL::Exception(oss.str().c_str());
1647 std::ostringstream oss; oss << "Cell # " << i << " of mesh \"" << mesh->getName() << "\" has no gauss location !";
1648 throw INTERP_KERNEL::Exception(oss.str().c_str());
1650 if(mesh->getTypeOfCell(i)!=_loc[dc[i]].getType())
1652 std::ostringstream oss; oss << "Types of mesh and gauss location mismatch for cell # " << i;
1653 throw INTERP_KERNEL::Exception(oss.str().c_str());
1656 int nbOfTuples=getNumberOfTuples(mesh);
1657 if(nbOfTuples!=da->getNumberOfTuples())
1659 std::ostringstream oss; oss << "Invalid number of tuples in the array : expecting " << nbOfTuples << " !";
1660 throw INTERP_KERNEL::Exception(oss.str().c_str());
1664 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationGauss::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
1667 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getMeasureField : mesh instance specified is NULL !");
1668 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=mesh->getMeasureField(isAbs);
1669 const double *volPtr=vol->getArray()->begin();
1670 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> ret=MEDCouplingFieldDouble::New(ON_GAUSS_PT);
1672 ret->setDiscretization(const_cast<MEDCouplingFieldDiscretizationGauss *>(this));
1673 if(!_discr_per_cell)
1674 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getMeasureField : no discr per cell array not defined ! spatial localization is incorrect !");
1675 _discr_per_cell->checkAllocated();
1676 if(_discr_per_cell->getNumberOfComponents()!=1)
1677 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::getMeasureField : no discr per cell array defined but with nb of components different from 1 !");
1678 if(_discr_per_cell->getNumberOfTuples()!=vol->getNumberOfTuples())
1679 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 !");
1680 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> offset=getOffsetArr(mesh);
1681 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arr=DataArrayDouble::New(); arr->alloc(getNumberOfTuples(mesh),1);
1683 double *arrPtr=arr->getPointer();
1684 const int *offsetPtr=offset->getConstPointer();
1685 int maxGaussLoc=(int)_loc.size();
1686 std::vector<int> locIds;
1687 std::vector<DataArrayInt *> ids=splitIntoSingleGaussDicrPerCellType(locIds);
1688 std::vector< MEDCouplingAutoRefCountObjectPtr<DataArrayInt> > ids2(ids.size()); std::copy(ids.begin(),ids.end(),ids2.begin());
1689 for(std::size_t i=0;i<locIds.size();i++)
1691 const DataArrayInt *curIds=ids[i];
1692 int locId=locIds[i];
1693 if(locId>=0 && locId<maxGaussLoc)
1695 const MEDCouplingGaussLocalization& loc=_loc[locId];
1696 int nbOfGaussPt=loc.getNumberOfGaussPt();
1697 INTERP_KERNEL::AutoPtr<double> weights=new double[nbOfGaussPt];
1698 double sum=std::accumulate(loc.getWeights().begin(),loc.getWeights().end(),0.);
1699 std::transform(loc.getWeights().begin(),loc.getWeights().end(),(double *)weights,std::bind2nd(std::multiplies<double>(),1./sum));
1700 for(const int *cellId=curIds->begin();cellId!=curIds->end();cellId++)
1701 for(int j=0;j<nbOfGaussPt;j++)
1702 arrPtr[offsetPtr[*cellId]+j]=weights[j]*volPtr[*cellId];
1706 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::getMeasureField : Presence of localization id " << locId << " in cell #" << curIds->getIJ(0,0) << " ! Must be in [0," << maxGaussLoc << ") !";
1707 throw INTERP_KERNEL::Exception(oss.str().c_str());
1710 ret->synchronizeTimeWithSupport();
1714 void MEDCouplingFieldDiscretizationGauss::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
1716 throw INTERP_KERNEL::Exception("Not implemented yet !");
1719 void MEDCouplingFieldDiscretizationGauss::getValueOnPos(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, int i, int j, int k, double *res) const
1721 throw INTERP_KERNEL::Exception("getValueOnPos(i,j,k) : Not applyable for Gauss points !");
1724 DataArrayDouble *MEDCouplingFieldDiscretizationGauss::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const
1726 throw INTERP_KERNEL::Exception("getValueOnMulti : Not implemented yet for gauss points !");
1729 MEDCouplingMesh *MEDCouplingFieldDiscretizationGauss::buildSubMeshData(const MEDCouplingMesh *mesh, const int *start, const int *end, DataArrayInt *&di) const
1732 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::buildSubMeshData : NULL input mesh !");
1733 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> diSafe=computeTupleIdsToSelectFromCellIds(mesh,start,end);
1734 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPart(start,end);
1740 * This method is strictly equivalent to MEDCouplingFieldDiscretizationGauss::buildSubMeshData except that it is optimized for input defined as a range of cell ids.
1742 * \param [out] beginOut Valid only if \a di is NULL
1743 * \param [out] endOut Valid only if \a di is NULL
1744 * \param [out] stepOut Valid only if \a di is NULL
1745 * \param [out] di is an array returned that specifies entity ids (nodes, cells, Gauss points... ) in array if no output range is foundable.
1747 * \sa MEDCouplingFieldDiscretizationGauss::buildSubMeshData
1749 MEDCouplingMesh *MEDCouplingFieldDiscretizationGauss::buildSubMeshDataRange(const MEDCouplingMesh *mesh, int beginCellIds, int endCellIds, int stepCellIds, int& beginOut, int& endOut, int& stepOut, DataArrayInt *&di) const
1751 if(stepCellIds!=1)//even for stepCellIds==-1 the output will not be a range
1752 return MEDCouplingFieldDiscretization::buildSubMeshDataRange(mesh,beginCellIds,endCellIds,stepCellIds,beginOut,endOut,stepOut,di);
1754 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::buildSubMeshDataRange : NULL input mesh !");
1755 if(!_discr_per_cell)
1756 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::buildSubMeshDataRange : no discretization array set !");
1757 di=0; beginOut=0; endOut=0; stepOut=stepCellIds;
1758 const char msg[]="MEDCouplingFieldDiscretizationGauss::buildSubMeshDataRange : cell #";
1759 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1760 const int *w=_discr_per_cell->begin();
1761 int nbMaxOfLocId=(int)_loc.size();
1762 for(int i=0;i<nbOfTuples;i++,w++)
1764 if(*w!=DFT_INVALID_LOCID_VALUE)
1766 if(*w>=0 && *w<nbMaxOfLocId)
1768 int delta=_loc[*w].getNumberOfGaussPt();
1776 { std::ostringstream oss; oss << msg << i << " has invalid id (" << *w << ") ! Should be in [0," << nbMaxOfLocId << ") !"; throw INTERP_KERNEL::Exception(oss.str().c_str()); }
1779 { std::ostringstream oss; oss << msg << i << " is detected as orphan !"; throw INTERP_KERNEL::Exception(oss.str().c_str()); }
1781 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPartRange(beginCellIds,endCellIds,stepCellIds);
1786 * This method returns a tuple ids selection from cell ids selection [start;end).
1787 * This method is called by MEDCouplingFieldDiscretizationGauss::buildSubMeshData to return parameter \b di.
1789 * \return a newly allocated array containing ids to select into the DataArrayDouble of the field.
1792 DataArrayInt *MEDCouplingFieldDiscretizationGauss::computeTupleIdsToSelectFromCellIds(const MEDCouplingMesh *mesh, const int *startCellIds, const int *endCellIds) const
1795 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::computeTupleIdsToSelectFromCellIds : null mesh !");
1796 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=buildNbOfGaussPointPerCellField();//check of _discr_per_cell not NULL pointer
1797 int nbOfCells=mesh->getNumberOfCells();
1798 if(_discr_per_cell->getNumberOfTuples()!=nbOfCells)
1799 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::computeTupleIdsToSelectFromCellIds : mismatch of nb of tuples of cell ids array and number of cells !");
1800 nbOfNodesPerCell->computeOffsets2();
1801 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> sel=DataArrayInt::New(); sel->useArray(startCellIds,false,CPP_DEALLOC,(int)std::distance(startCellIds,endCellIds),1);
1802 return sel->buildExplicitArrByRanges(nbOfNodesPerCell);
1806 * No implementation needed !
1808 void MEDCouplingFieldDiscretizationGauss::renumberValuesOnNodes(double , const int *, int newNbOfNodes, DataArrayDouble *) const
1812 void MEDCouplingFieldDiscretizationGauss::renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const
1814 throw INTERP_KERNEL::Exception("Not implemented yet !");
1817 void MEDCouplingFieldDiscretizationGauss::renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const
1819 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 !");
1822 void MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnType(const MEDCouplingMesh *mesh, INTERP_KERNEL::NormalizedCellType type, const std::vector<double>& refCoo,
1823 const std::vector<double>& gsCoo, const std::vector<double>& wg)
1826 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnType : NULL input mesh !");
1827 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
1828 if((int)cm.getDimension()!=mesh->getMeshDimension())
1830 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnType : mismatch of dimensions ! MeshDim==" << mesh->getMeshDimension();
1831 oss << " whereas Type '" << cm.getRepr() << "' has dimension " << cm.getDimension() << " !";
1832 throw INTERP_KERNEL::Exception(oss.str().c_str());
1834 buildDiscrPerCellIfNecessary(mesh);
1835 int id=(int)_loc.size();
1836 MEDCouplingGaussLocalization elt(type,refCoo,gsCoo,wg);
1837 _loc.push_back(elt);
1838 int *ptr=_discr_per_cell->getPointer();
1839 int nbCells=mesh->getNumberOfCells();
1840 for(int i=0;i<nbCells;i++)
1841 if(mesh->getTypeOfCell(i)==type)
1843 zipGaussLocalizations();
1846 void MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnCells(const MEDCouplingMesh *mesh, const int *begin, const int *end, const std::vector<double>& refCoo,
1847 const std::vector<double>& gsCoo, const std::vector<double>& wg)
1850 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnCells : NULL input mesh !");
1851 buildDiscrPerCellIfNecessary(mesh);
1852 if(std::distance(begin,end)<1)
1853 throw INTERP_KERNEL::Exception("Size of [begin,end) must be equal or greater than 1 !");
1854 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(*begin);
1855 MEDCouplingGaussLocalization elt(type,refCoo,gsCoo,wg);
1856 int id=(int)_loc.size();
1857 int *ptr=_discr_per_cell->getPointer();
1858 for(const int *w=begin+1;w!=end;w++)
1860 if(mesh->getTypeOfCell(*w)!=type)
1862 std::ostringstream oss; oss << "The cell with id " << *w << " has been detected to be incompatible in the [begin,end) array specified !";
1863 throw INTERP_KERNEL::Exception(oss.str().c_str());
1867 for(const int *w2=begin;w2!=end;w2++)
1870 _loc.push_back(elt);
1871 zipGaussLocalizations();
1874 void MEDCouplingFieldDiscretizationGauss::clearGaussLocalizations()
1878 _discr_per_cell->decrRef();
1884 void MEDCouplingFieldDiscretizationGauss::setGaussLocalization(int locId, const MEDCouplingGaussLocalization& loc)
1887 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::setGaussLocalization : localization id has to be >=0 !");
1888 int sz=(int)_loc.size();
1889 MEDCouplingGaussLocalization gLoc(INTERP_KERNEL::NORM_ERROR);
1891 _loc.resize(locId+1,gLoc);
1895 void MEDCouplingFieldDiscretizationGauss::resizeLocalizationVector(int newSz)
1898 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::resizeLocalizationVector : new size has to be >=0 !");
1899 MEDCouplingGaussLocalization gLoc(INTERP_KERNEL::NORM_ERROR);
1900 _loc.resize(newSz,gLoc);
1903 MEDCouplingGaussLocalization& MEDCouplingFieldDiscretizationGauss::getGaussLocalization(int locId)
1905 checkLocalizationId(locId);
1909 int MEDCouplingFieldDiscretizationGauss::getNbOfGaussLocalization() const
1911 return (int)_loc.size();
1914 int MEDCouplingFieldDiscretizationGauss::getGaussLocalizationIdOfOneCell(int cellId) const
1916 if(!_discr_per_cell)
1917 throw INTERP_KERNEL::Exception("No Gauss localization still set !");
1918 int locId=_discr_per_cell->begin()[cellId];
1920 throw INTERP_KERNEL::Exception("No Gauss localization set for the specified cell !");
1924 int MEDCouplingFieldDiscretizationGauss::getGaussLocalizationIdOfOneType(INTERP_KERNEL::NormalizedCellType type) const
1926 std::set<int> ret=getGaussLocalizationIdsOfOneType(type);
1928 throw INTERP_KERNEL::Exception("No gauss discretization found for the specified type !");
1930 throw INTERP_KERNEL::Exception("Several gauss discretizations have been found for the specified type !");
1931 return *ret.begin();
1934 std::set<int> MEDCouplingFieldDiscretizationGauss::getGaussLocalizationIdsOfOneType(INTERP_KERNEL::NormalizedCellType type) const
1936 if(!_discr_per_cell)
1937 throw INTERP_KERNEL::Exception("No Gauss localization still set !");
1940 for(std::vector<MEDCouplingGaussLocalization>::const_iterator iter=_loc.begin();iter!=_loc.end();iter++,id++)
1941 if((*iter).getType()==type)
1946 void MEDCouplingFieldDiscretizationGauss::getCellIdsHavingGaussLocalization(int locId, std::vector<int>& cellIds) const
1948 if(locId<0 || locId>=(int)_loc.size())
1949 throw INTERP_KERNEL::Exception("Invalid locId given : must be in range [0:getNbOfGaussLocalization()) !");
1950 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1951 const int *ptr=_discr_per_cell->getConstPointer();
1952 for(int i=0;i<nbOfTuples;i++)
1954 cellIds.push_back(i);
1957 const MEDCouplingGaussLocalization& MEDCouplingFieldDiscretizationGauss::getGaussLocalization(int locId) const
1959 checkLocalizationId(locId);
1963 void MEDCouplingFieldDiscretizationGauss::checkLocalizationId(int locId) const
1965 if(locId<0 || locId>=(int)_loc.size())
1966 throw INTERP_KERNEL::Exception("Invalid locId given : must be in range [0:getNbOfGaussLocalization()) !");
1969 int MEDCouplingFieldDiscretizationGauss::getOffsetOfCell(int cellId) const
1972 const int *start=_discr_per_cell->getConstPointer();
1973 for(const int *w=start;w!=start+cellId;w++)
1974 ret+=_loc[*w].getNumberOfGaussPt();
1979 * This method do the assumption that there is no orphan cell. If there is an exception is thrown.
1980 * This method makes the assumption too that '_discr_per_cell' is defined. If not an exception is thrown.
1981 * This method returns a newly created array with number of tuples equals to '_discr_per_cell->getNumberOfTuples' and number of components equal to 1.
1982 * The i_th tuple in returned array is the number of gauss point if the corresponding cell.
1984 DataArrayInt *MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField() const
1986 if(!_discr_per_cell)
1987 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField : no discretization array set !");
1988 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
1989 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret=DataArrayInt::New();
1990 const int *w=_discr_per_cell->begin();
1991 ret->alloc(nbOfTuples,1);
1992 int *valsToFill=ret->getPointer();
1993 int nbMaxOfLocId=(int)_loc.size();
1994 for(int i=0;i<nbOfTuples;i++,w++)
1995 if(*w!=DFT_INVALID_LOCID_VALUE)
1997 if(*w>=0 && *w<nbMaxOfLocId)
1998 valsToFill[i]=_loc[*w].getNumberOfGaussPt();
2001 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField : cell #" << i << " has invalid id (" << *w << ") ! Should be in [0," << nbMaxOfLocId << ") !";
2002 throw INTERP_KERNEL::Exception(oss.str().c_str());
2007 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGauss::buildNbOfGaussPointPerCellField : cell #" << i << " is detected as orphan !";
2008 throw INTERP_KERNEL::Exception(oss.str().c_str());
2013 void MEDCouplingFieldDiscretizationGauss::reprQuickOverview(std::ostream& stream) const
2015 stream << "Gauss points spatial discretization.";
2019 * This method makes the assumption that _discr_per_cell is set.
2020 * This method reduces as much as possible number size of _loc.
2021 * This method is useful when several set on same cells has been done and that some Gauss Localization are no more used.
2023 void MEDCouplingFieldDiscretizationGauss::zipGaussLocalizations()
2025 const int *start=_discr_per_cell->begin();
2026 int nbOfTuples=_discr_per_cell->getNumberOfTuples();
2027 INTERP_KERNEL::AutoPtr<int> tmp=new int[_loc.size()];
2028 std::fill((int *)tmp,(int *)tmp+_loc.size(),-2);
2029 for(const int *w=start;w!=start+nbOfTuples;w++)
2033 for(int i=0;i<(int)_loc.size();i++)
2036 if(fid==(int)_loc.size())
2039 int *start2=_discr_per_cell->getPointer();
2040 for(int *w2=start2;w2!=start2+nbOfTuples;w2++)
2043 std::vector<MEDCouplingGaussLocalization> tmpLoc;
2044 for(int i=0;i<(int)_loc.size();i++)
2046 tmpLoc.push_back(_loc[i]);
2050 MEDCouplingFieldDiscretizationGaussNE::MEDCouplingFieldDiscretizationGaussNE()
2054 TypeOfField MEDCouplingFieldDiscretizationGaussNE::getEnum() const
2060 * This method is simply called by MEDCouplingFieldDiscretization::deepCpy. It performs the deep copy of \a this.
2062 * \sa MEDCouplingFieldDiscretization::deepCpy.
2064 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationGaussNE::clone() const
2066 return new MEDCouplingFieldDiscretizationGaussNE(*this);
2069 std::string MEDCouplingFieldDiscretizationGaussNE::getStringRepr() const
2071 return std::string(REPR);
2074 const char *MEDCouplingFieldDiscretizationGaussNE::getRepr() const
2079 bool MEDCouplingFieldDiscretizationGaussNE::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
2083 reason="other spatial discretization is NULL, and this spatial discretization (GaussNE) is defined.";
2086 const MEDCouplingFieldDiscretizationGaussNE *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationGaussNE *>(other);
2089 reason="Spatial discrtization of this is ON_GAUSS_NE, which is not the case of other.";
2094 * This method returns the number of tuples regarding exclusively the input code \b without \b using \b a \b mesh \b in \b input.
2095 * The input code coherency is also checked regarding spatial discretization of \a this.
2096 * If an incoherency is detected, an exception will be thrown. If the input code is coherent, the number of tuples expected is returned.
2097 * 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).
2099 int MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode(const std::vector<int>& code, const std::vector<const DataArrayInt *>& idsPerType) const
2101 if(code.size()%3!=0)
2102 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode : invalid input code !");
2103 int nbOfSplit=(int)idsPerType.size();
2104 int nbOfTypes=(int)code.size()/3;
2106 for(int i=0;i<nbOfTypes;i++)
2108 const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel((INTERP_KERNEL::NormalizedCellType)code[3*i]));
2111 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 !";
2112 throw INTERP_KERNEL::Exception(oss.str().c_str());
2114 int nbOfEltInChunk=code[3*i+1];
2115 if(nbOfEltInChunk<0)
2116 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode : invalid input code ! presence of negative value in a type !");
2117 int pos=code[3*i+2];
2120 if(pos<0 || pos>=nbOfSplit)
2122 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode : input code points to pos " << pos << " in typeid " << i << " ! Should be in [0," << nbOfSplit << ") !";
2123 throw INTERP_KERNEL::Exception(oss.str().c_str());
2125 const DataArrayInt *ids(idsPerType[pos]);
2126 if(!ids || !ids->isAllocated() || ids->getNumberOfComponents()!=1 || ids->getNumberOfTuples()!=nbOfEltInChunk || ids->getMinValueInArray()<0)
2128 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuplesExpectedRegardingCode : input pfl chunck at pos " << pos << " should have " << i << " tuples and one component and with ids all >=0 !";
2129 throw INTERP_KERNEL::Exception(oss.str().c_str());
2132 ret+=nbOfEltInChunk*(int)cm.getNumberOfNodes();
2137 int MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuples(const MEDCouplingMesh *mesh) const
2140 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getNumberOfTuples : NULL input mesh !");
2142 int nbOfCells=mesh->getNumberOfCells();
2143 for(int i=0;i<nbOfCells;i++)
2145 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
2146 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2148 throw INTERP_KERNEL::Exception("Not implemented yet Gauss node on elements for polygons and polyedrons !");
2149 ret+=cm.getNumberOfNodes();
2154 int MEDCouplingFieldDiscretizationGaussNE::getNumberOfMeshPlaces(const MEDCouplingMesh *mesh) const
2157 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getNumberOfMeshPlaces : NULL input mesh !");
2158 return mesh->getNumberOfCells();
2161 DataArrayInt *MEDCouplingFieldDiscretizationGaussNE::getOffsetArr(const MEDCouplingMesh *mesh) const
2164 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getOffsetArr : NULL input mesh !");
2165 int nbOfTuples=mesh->getNumberOfCells();
2166 DataArrayInt *ret=DataArrayInt::New();
2167 ret->alloc(nbOfTuples+1,1);
2168 int *retPtr=ret->getPointer();
2170 for(int i=0;i<nbOfTuples;i++)
2172 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
2173 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2175 throw INTERP_KERNEL::Exception("Not implemented yet Gauss node on elements for polygons and polyedrons !");
2176 retPtr[i+1]=retPtr[i]+cm.getNumberOfNodes();
2181 void MEDCouplingFieldDiscretizationGaussNE::renumberArraysForCell(const MEDCouplingMesh *mesh, const std::vector<DataArray *>& arrays,
2182 const int *old2NewBg, bool check)
2185 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::renumberArraysForCell : NULL input mesh !");
2186 const int *array=old2NewBg;
2188 array=DataArrayInt::CheckAndPreparePermutation(old2NewBg,old2NewBg+mesh->getNumberOfCells());
2189 int nbOfCells=mesh->getNumberOfCells();
2190 int nbOfTuples=getNumberOfTuples(mesh);
2191 int *array2=new int[nbOfTuples];//stores the final conversion array old2New to give to arrays in renumberInPlace.
2192 int *array3=new int[nbOfCells];//store for each cell in after renumbering the offset needed by each cell in new numbering.
2194 for(int i=1;i<nbOfCells;i++)
2196 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell((int)std::distance(array,std::find(array,array+nbOfCells,i-1)));
2197 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2198 array3[i]=array3[i-1]+cm.getNumberOfNodes();
2201 for(int i=0;i<nbOfCells;i++)
2203 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
2204 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2205 for(int k=0;k<(int)cm.getNumberOfNodes();k++,j++)
2206 array2[j]=array3[array[i]]+k;
2209 for(std::vector<DataArray *>::const_iterator it=arrays.begin();it!=arrays.end();it++)
2211 (*it)->renumberInPlace(array2);
2214 free(const_cast<int *>(array));
2217 DataArrayDouble *MEDCouplingFieldDiscretizationGaussNE::getLocalizationOfDiscValues(const MEDCouplingMesh *mesh) const
2220 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getLocalizationOfDiscValues : NULL input mesh !");
2221 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
2222 MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> umesh=mesh->buildUnstructured();//in general do nothing
2223 int nbOfTuples=getNumberOfTuples(umesh);
2224 int spaceDim=mesh->getSpaceDimension();
2225 ret->alloc(nbOfTuples,spaceDim);
2226 const double *coords=umesh->getCoords()->begin();
2227 const int *connI=umesh->getNodalConnectivityIndex()->getConstPointer();
2228 const int *conn=umesh->getNodalConnectivity()->getConstPointer();
2229 int nbCells=umesh->getNumberOfCells();
2230 double *retPtr=ret->getPointer();
2231 for(int i=0;i<nbCells;i++,connI++)
2232 for(const int *w=conn+connI[0]+1;w!=conn+connI[1];w++)
2234 retPtr=std::copy(coords+(*w)*spaceDim,coords+((*w)+1)*spaceDim,retPtr);
2239 * Reimplemented from MEDCouplingFieldDiscretization::integral for performance reason. The default implementation is valid too for GAUSS_NE spatial discretization.
2241 void MEDCouplingFieldDiscretizationGaussNE::integral(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, bool isWAbs, double *res) const
2244 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::integral : input mesh or array is null !");
2245 int nbOfCompo=arr->getNumberOfComponents();
2246 std::fill(res,res+nbOfCompo,0.);
2248 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=mesh->getMeasureField(isWAbs);
2249 std::set<INTERP_KERNEL::NormalizedCellType> types=mesh->getAllGeoTypes();
2250 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=mesh->computeNbOfNodesPerCell();
2251 nbOfNodesPerCell->computeOffsets2();
2252 const double *arrPtr=arr->begin(),*volPtr=vol->getArray()->begin();
2253 for(std::set<INTERP_KERNEL::NormalizedCellType>::const_iterator it=types.begin();it!=types.end();it++)
2255 std::size_t wArrSz=-1;
2256 const double *wArr=GetWeightArrayFromGeometricType(*it,wArrSz);
2257 INTERP_KERNEL::AutoPtr<double> wArr2=new double[wArrSz];
2258 double sum=std::accumulate(wArr,wArr+wArrSz,0.);
2259 std::transform(wArr,wArr+wArrSz,(double *)wArr2,std::bind2nd(std::multiplies<double>(),1./sum));
2260 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ids=mesh->giveCellsWithType(*it);
2261 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ids2=ids->buildExplicitArrByRanges(nbOfNodesPerCell);
2262 const int *ptIds2=ids2->begin(),*ptIds=ids->begin();
2263 int nbOfCellsWithCurGeoType=ids->getNumberOfTuples();
2264 for(int i=0;i<nbOfCellsWithCurGeoType;i++,ptIds++,ptIds2+=wArrSz)
2266 for(int k=0;k<nbOfCompo;k++)
2269 for(std::size_t j=0;j<wArrSz;j++)
2270 tmp+=arrPtr[nbOfCompo*ptIds2[j]+k]*wArr2[j];
2271 res[k]+=tmp*volPtr[*ptIds];
2277 const double *MEDCouplingFieldDiscretizationGaussNE::GetWeightArrayFromGeometricType(INTERP_KERNEL::NormalizedCellType geoType, std::size_t& lgth)
2281 case INTERP_KERNEL::NORM_POINT1:
2282 lgth=(int)sizeof(FGP_POINT1)/sizeof(double);
2284 case INTERP_KERNEL::NORM_SEG2:
2285 lgth=(int)sizeof(FGP_SEG2)/sizeof(double);
2287 case INTERP_KERNEL::NORM_SEG3:
2288 lgth=(int)sizeof(FGP_SEG3)/sizeof(double);
2290 case INTERP_KERNEL::NORM_SEG4:
2291 lgth=(int)sizeof(FGP_SEG4)/sizeof(double);
2293 case INTERP_KERNEL::NORM_TRI3:
2294 lgth=(int)sizeof(FGP_TRI3)/sizeof(double);
2296 case INTERP_KERNEL::NORM_TRI6:
2297 lgth=(int)sizeof(FGP_TRI6)/sizeof(double);
2299 case INTERP_KERNEL::NORM_TRI7:
2300 lgth=(int)sizeof(FGP_TRI7)/sizeof(double);
2302 case INTERP_KERNEL::NORM_QUAD4:
2303 lgth=(int)sizeof(FGP_QUAD4)/sizeof(double);
2305 case INTERP_KERNEL::NORM_QUAD8:
2306 lgth=(int)sizeof(FGP_QUAD8)/sizeof(double);
2308 case INTERP_KERNEL::NORM_QUAD9:
2309 lgth=(int)sizeof(FGP_QUAD9)/sizeof(double);
2311 case INTERP_KERNEL::NORM_TETRA4:
2312 lgth=(int)sizeof(FGP_TETRA4)/sizeof(double);
2314 case INTERP_KERNEL::NORM_TETRA10:
2315 lgth=(int)sizeof(FGP_TETRA10)/sizeof(double);
2317 case INTERP_KERNEL::NORM_PENTA6:
2318 lgth=(int)sizeof(FGP_PENTA6)/sizeof(double);
2320 case INTERP_KERNEL::NORM_PENTA15:
2321 lgth=(int)sizeof(FGP_PENTA15)/sizeof(double);
2323 case INTERP_KERNEL::NORM_HEXA8:
2324 lgth=(int)sizeof(FGP_HEXA8)/sizeof(double);
2326 case INTERP_KERNEL::NORM_HEXA20:
2327 lgth=(int)sizeof(FGP_HEXA20)/sizeof(double);
2329 case INTERP_KERNEL::NORM_HEXA27:
2330 lgth=(int)sizeof(FGP_HEXA27)/sizeof(double);
2332 case INTERP_KERNEL::NORM_PYRA5:
2333 lgth=(int)sizeof(FGP_PYRA5)/sizeof(double);
2335 case INTERP_KERNEL::NORM_PYRA13:
2336 lgth=(int)sizeof(FGP_PYRA13)/sizeof(double);
2339 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::GetWeightArrayFromGeometricType : only SEG[2,3,4], TRI[3,6,7], QUAD[4,9], TETRA[4,10], PENTA[6,15], HEXA[8,20,27], PYRA[5,13] supported !");
2343 const double *MEDCouplingFieldDiscretizationGaussNE::GetRefCoordsFromGeometricType(INTERP_KERNEL::NormalizedCellType geoType, std::size_t& lgth)
2347 case INTERP_KERNEL::NORM_POINT1:
2350 case INTERP_KERNEL::NORM_SEG2:
2351 lgth=(int)sizeof(REF_SEG2)/sizeof(double);
2353 case INTERP_KERNEL::NORM_SEG3:
2354 lgth=(int)sizeof(REF_SEG3)/sizeof(double);
2356 case INTERP_KERNEL::NORM_SEG4:
2357 lgth=(int)sizeof(REF_SEG4)/sizeof(double);
2359 case INTERP_KERNEL::NORM_TRI3:
2360 lgth=(int)sizeof(REF_TRI3)/sizeof(double);
2362 case INTERP_KERNEL::NORM_TRI6:
2363 lgth=(int)sizeof(REF_TRI6)/sizeof(double);
2365 case INTERP_KERNEL::NORM_TRI7:
2366 lgth=(int)sizeof(REF_TRI7)/sizeof(double);
2368 case INTERP_KERNEL::NORM_QUAD4:
2369 lgth=(int)sizeof(REF_QUAD4)/sizeof(double);
2371 case INTERP_KERNEL::NORM_QUAD8:
2372 lgth=(int)sizeof(REF_QUAD8)/sizeof(double);
2374 case INTERP_KERNEL::NORM_QUAD9:
2375 lgth=(int)sizeof(REF_QUAD9)/sizeof(double);
2377 case INTERP_KERNEL::NORM_TETRA4:
2378 lgth=(int)sizeof(REF_TETRA4)/sizeof(double);
2380 case INTERP_KERNEL::NORM_TETRA10:
2381 lgth=(int)sizeof(REF_TETRA10)/sizeof(double);
2383 case INTERP_KERNEL::NORM_PENTA6:
2384 lgth=(int)sizeof(REF_PENTA6)/sizeof(double);
2386 case INTERP_KERNEL::NORM_PENTA15:
2387 lgth=(int)sizeof(REF_PENTA15)/sizeof(double);
2389 case INTERP_KERNEL::NORM_HEXA8:
2390 lgth=(int)sizeof(REF_HEXA8)/sizeof(double);
2392 case INTERP_KERNEL::NORM_HEXA20:
2393 lgth=(int)sizeof(REF_HEXA20)/sizeof(double);
2395 case INTERP_KERNEL::NORM_HEXA27:
2396 lgth=(int)sizeof(REF_HEXA27)/sizeof(double);
2398 case INTERP_KERNEL::NORM_PYRA5:
2399 lgth=(int)sizeof(REF_PYRA5)/sizeof(double);
2401 case INTERP_KERNEL::NORM_PYRA13:
2402 lgth=(int)sizeof(REF_PYRA13)/sizeof(double);
2405 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 !");
2409 const double *MEDCouplingFieldDiscretizationGaussNE::GetLocsFromGeometricType(INTERP_KERNEL::NormalizedCellType geoType, std::size_t& lgth)
2413 case INTERP_KERNEL::NORM_POINT1:
2418 case INTERP_KERNEL::NORM_SEG2:
2420 lgth=(int)sizeof(LOC_SEG2)/sizeof(double);
2423 case INTERP_KERNEL::NORM_SEG3:
2425 lgth=(int)sizeof(LOC_SEG3)/sizeof(double);
2428 case INTERP_KERNEL::NORM_SEG4:
2430 lgth=(int)sizeof(LOC_SEG4)/sizeof(double);
2433 case INTERP_KERNEL::NORM_TRI3:
2435 lgth=(int)sizeof(LOC_TRI3)/sizeof(double);
2438 case INTERP_KERNEL::NORM_TRI6:
2440 lgth=(int)sizeof(LOC_TRI6)/sizeof(double);
2443 case INTERP_KERNEL::NORM_TRI7:
2445 lgth=(int)sizeof(LOC_TRI7)/sizeof(double);
2448 case INTERP_KERNEL::NORM_QUAD4:
2450 lgth=(int)sizeof(LOC_QUAD4)/sizeof(double);
2453 case INTERP_KERNEL::NORM_QUAD8:
2455 lgth=(int)sizeof(LOC_QUAD8)/sizeof(double);
2458 case INTERP_KERNEL::NORM_QUAD9:
2460 lgth=(int)sizeof(LOC_QUAD9)/sizeof(double);
2463 case INTERP_KERNEL::NORM_TETRA4:
2465 lgth=(int)sizeof(LOC_TETRA4)/sizeof(double);
2468 case INTERP_KERNEL::NORM_TETRA10:
2470 lgth=(int)sizeof(LOC_TETRA10)/sizeof(double);
2473 case INTERP_KERNEL::NORM_PENTA6:
2475 lgth=(int)sizeof(LOC_PENTA6)/sizeof(double);
2478 case INTERP_KERNEL::NORM_PENTA15:
2480 lgth=(int)sizeof(LOC_PENTA15)/sizeof(double);
2483 case INTERP_KERNEL::NORM_HEXA8:
2485 lgth=(int)sizeof(LOC_HEXA8)/sizeof(double);
2488 case INTERP_KERNEL::NORM_HEXA20:
2490 lgth=(int)sizeof(LOC_HEXA20)/sizeof(double);
2493 case INTERP_KERNEL::NORM_HEXA27:
2495 lgth=(int)sizeof(LOC_HEXA27)/sizeof(double);
2498 case INTERP_KERNEL::NORM_PYRA5:
2500 lgth=(int)sizeof(LOC_PYRA5)/sizeof(double);
2503 case INTERP_KERNEL::NORM_PYRA13:
2505 lgth=(int)sizeof(LOC_PYRA13)/sizeof(double);
2509 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 !");
2513 void MEDCouplingFieldDiscretizationGaussNE::computeMeshRestrictionFromTupleIds(const MEDCouplingMesh *mesh, const int *tupleIdsBg, const int *tupleIdsEnd,
2514 DataArrayInt *&cellRestriction, DataArrayInt *&trueTupleRestriction) const
2517 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::computeMeshRestrictionFromTupleIds : NULL input mesh !");
2518 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> tmp=DataArrayInt::New(); tmp->alloc((int)std::distance(tupleIdsBg,tupleIdsEnd),1);
2519 std::copy(tupleIdsBg,tupleIdsEnd,tmp->getPointer());
2521 tmp=tmp->buildUnique();
2522 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=mesh->computeNbOfNodesPerCell();
2523 nbOfNodesPerCell->computeOffsets2();
2524 nbOfNodesPerCell->searchRangesInListOfIds(tmp,cellRestriction,trueTupleRestriction);
2527 void MEDCouplingFieldDiscretizationGaussNE::checkCompatibilityWithNature(NatureOfField nat) const
2531 double MEDCouplingFieldDiscretizationGaussNE::getIJK(const MEDCouplingMesh *mesh, const DataArrayDouble *da, int cellId, int nodeIdInCell, int compoId) const
2534 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getIJK : NULL input mesh !");
2536 for(int i=0;i<cellId;i++)
2538 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
2539 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2540 offset+=cm.getNumberOfNodes();
2542 return da->getIJ(offset+nodeIdInCell,compoId);
2545 void MEDCouplingFieldDiscretizationGaussNE::checkCoherencyBetween(const MEDCouplingMesh *mesh, const DataArray *da) const
2547 int nbOfTuples=getNumberOfTuples(mesh);
2548 if(nbOfTuples!=da->getNumberOfTuples())
2550 std::ostringstream oss; oss << "Invalid number of tuples in the array : expecting " << nbOfTuples << " !";
2551 throw INTERP_KERNEL::Exception(oss.str().c_str());
2555 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationGaussNE::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
2558 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::getMeasureField : mesh instance specified is NULL !");
2559 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> vol=mesh->getMeasureField(isAbs);
2560 const double *volPtr=vol->getArray()->begin();
2561 MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> ret=MEDCouplingFieldDouble::New(ON_GAUSS_NE);
2564 std::set<INTERP_KERNEL::NormalizedCellType> types=mesh->getAllGeoTypes();
2565 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=mesh->computeNbOfNodesPerCell();
2566 int nbTuples=nbOfNodesPerCell->accumulate(0);
2567 nbOfNodesPerCell->computeOffsets2();
2568 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arr=DataArrayDouble::New(); arr->alloc(nbTuples,1);
2570 double *arrPtr=arr->getPointer();
2571 for(std::set<INTERP_KERNEL::NormalizedCellType>::const_iterator it=types.begin();it!=types.end();it++)
2573 std::size_t wArrSz=-1;
2574 const double *wArr=GetWeightArrayFromGeometricType(*it,wArrSz);
2575 INTERP_KERNEL::AutoPtr<double> wArr2=new double[wArrSz];
2576 double sum=std::accumulate(wArr,wArr+wArrSz,0.);
2577 std::transform(wArr,wArr+wArrSz,(double *)wArr2,std::bind2nd(std::multiplies<double>(),1./sum));
2578 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ids=mesh->giveCellsWithType(*it);
2579 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ids2=ids->buildExplicitArrByRanges(nbOfNodesPerCell);
2580 const int *ptIds2=ids2->begin(),*ptIds=ids->begin();
2581 int nbOfCellsWithCurGeoType=ids->getNumberOfTuples();
2582 for(int i=0;i<nbOfCellsWithCurGeoType;i++,ptIds++)
2583 for(std::size_t j=0;j<wArrSz;j++,ptIds2++)
2584 arrPtr[*ptIds2]=wArr2[j]*volPtr[*ptIds];
2586 ret->synchronizeTimeWithSupport();
2590 void MEDCouplingFieldDiscretizationGaussNE::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
2592 throw INTERP_KERNEL::Exception("Not implemented yet !");
2595 void MEDCouplingFieldDiscretizationGaussNE::getValueOnPos(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, int i, int j, int k, double *res) const
2597 throw INTERP_KERNEL::Exception("getValueOnPos(i,j,k) : Not applyable for Gauss points !");
2600 DataArrayDouble *MEDCouplingFieldDiscretizationGaussNE::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const
2602 throw INTERP_KERNEL::Exception("getValueOnMulti : Not implemented for Gauss NE !");
2605 MEDCouplingMesh *MEDCouplingFieldDiscretizationGaussNE::buildSubMeshData(const MEDCouplingMesh *mesh, const int *start, const int *end, DataArrayInt *&di) const
2608 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::buildSubMeshData : NULL input mesh !");
2609 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> diSafe=computeTupleIdsToSelectFromCellIds(mesh,start,end);
2610 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPart(start,end);
2616 * This method is strictly equivalent to MEDCouplingFieldDiscretizationGauss::buildSubMeshData except that it is optimized for input defined as a range of cell ids.
2618 * \param [out] beginOut Valid only if \a di is NULL
2619 * \param [out] endOut Valid only if \a di is NULL
2620 * \param [out] stepOut Valid only if \a di is NULL
2621 * \param [out] di is an array returned that specifies entity ids (nodes, cells, Gauss points... ) in array if no output range is foundable.
2623 * \sa MEDCouplingFieldDiscretizationGauss::buildSubMeshData
2625 MEDCouplingMesh *MEDCouplingFieldDiscretizationGaussNE::buildSubMeshDataRange(const MEDCouplingMesh *mesh, int beginCellIds, int endCellIds, int stepCellIds, int& beginOut, int& endOut, int& stepOut, DataArrayInt *&di) const
2627 if(stepCellIds!=1)//even for stepCellIds==-1 the output will not be a range
2628 return MEDCouplingFieldDiscretization::buildSubMeshDataRange(mesh,beginCellIds,endCellIds,stepCellIds,beginOut,endOut,stepOut,di);
2630 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::buildSubMeshDataRange : NULL input mesh !");
2631 int nbOfCells=mesh->getNumberOfCells();
2632 di=0; beginOut=0; endOut=0; stepOut=stepCellIds;
2633 const char msg[]="MEDCouplingFieldDiscretizationGaussNE::buildSubMeshDataRange : cell #";
2634 for(int i=0;i<nbOfCells;i++)
2636 INTERP_KERNEL::NormalizedCellType type=mesh->getTypeOfCell(i);
2637 const INTERP_KERNEL::CellModel& cm=INTERP_KERNEL::CellModel::GetCellModel(type);
2639 { std::ostringstream oss; oss << msg << i << " presence of dynamic cell (polygons and polyedrons) ! Not implemented !"; throw INTERP_KERNEL::Exception(oss.str().c_str()); }
2640 int delta=cm.getNumberOfNodes();
2647 MEDCouplingAutoRefCountObjectPtr<MEDCouplingMesh> ret=mesh->buildPartRange(beginCellIds,endCellIds,stepCellIds);
2653 * This method returns a tuple ids selection from cell ids selection [start;end).
2654 * This method is called by MEDCouplingFieldDiscretizationGaussNE::buildSubMeshData to return parameter \b di.
2656 * \return a newly allocated array containing ids to select into the DataArrayDouble of the field.
2659 DataArrayInt *MEDCouplingFieldDiscretizationGaussNE::computeTupleIdsToSelectFromCellIds(const MEDCouplingMesh *mesh, const int *startCellIds, const int *endCellIds) const
2662 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGaussNE::computeTupleIdsToSelectFromCellIds : null mesh !");
2663 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfNodesPerCell=mesh->computeNbOfNodesPerCell();
2664 nbOfNodesPerCell->computeOffsets2();
2665 MEDCouplingAutoRefCountObjectPtr<DataArrayInt> sel=DataArrayInt::New(); sel->useArray(startCellIds,false,CPP_DEALLOC,(int)std::distance(startCellIds,endCellIds),1);
2666 return sel->buildExplicitArrByRanges(nbOfNodesPerCell);
2670 * No implementation needed !
2672 void MEDCouplingFieldDiscretizationGaussNE::renumberValuesOnNodes(double , const int *, int newNbOfNodes, DataArrayDouble *) const
2676 void MEDCouplingFieldDiscretizationGaussNE::renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const
2678 throw INTERP_KERNEL::Exception("Not implemented yet !");
2681 void MEDCouplingFieldDiscretizationGaussNE::renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const
2683 throw INTERP_KERNEL::Exception("Not implemented yet !");
2686 void MEDCouplingFieldDiscretizationGaussNE::reprQuickOverview(std::ostream& stream) const
2688 stream << "Gauss points on nodes per element spatial discretization.";
2691 MEDCouplingFieldDiscretizationGaussNE::MEDCouplingFieldDiscretizationGaussNE(const MEDCouplingFieldDiscretizationGaussNE& other):MEDCouplingFieldDiscretization(other)
2695 TypeOfField MEDCouplingFieldDiscretizationKriging::getEnum() const
2700 const char *MEDCouplingFieldDiscretizationKriging::getRepr() const
2706 * This method is simply called by MEDCouplingFieldDiscretization::deepCpy. It performs the deep copy of \a this.
2708 * \sa MEDCouplingFieldDiscretization::deepCpy.
2710 MEDCouplingFieldDiscretization *MEDCouplingFieldDiscretizationKriging::clone() const
2712 return new MEDCouplingFieldDiscretizationKriging;
2715 std::string MEDCouplingFieldDiscretizationKriging::getStringRepr() const
2717 return std::string(REPR);
2720 void MEDCouplingFieldDiscretizationKriging::checkCompatibilityWithNature(NatureOfField nat) const
2722 if(nat!=ConservativeVolumic)
2723 throw INTERP_KERNEL::Exception("Invalid nature for Kriging field : expected ConservativeVolumic !");
2726 bool MEDCouplingFieldDiscretizationKriging::isEqualIfNotWhy(const MEDCouplingFieldDiscretization *other, double eps, std::string& reason) const
2730 reason="other spatial discretization is NULL, and this spatial discretization (Kriginig) is defined.";
2733 const MEDCouplingFieldDiscretizationKriging *otherC=dynamic_cast<const MEDCouplingFieldDiscretizationKriging *>(other);
2736 reason="Spatial discrtization of this is ON_NODES_KR, which is not the case of other.";
2740 MEDCouplingFieldDouble *MEDCouplingFieldDiscretizationKriging::getMeasureField(const MEDCouplingMesh *mesh, bool isAbs) const
2743 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationKriging::getMeasureField : mesh instance specified is NULL !");
2744 throw INTERP_KERNEL::Exception("getMeasureField on FieldDiscretizationKriging : not implemented yet !");
2747 void MEDCouplingFieldDiscretizationKriging::getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const
2749 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> res2=MEDCouplingFieldDiscretizationKriging::getValueOnMulti(arr,mesh,loc,1);
2750 std::copy(res2->begin(),res2->end(),res);
2753 DataArrayDouble *MEDCouplingFieldDiscretizationKriging::getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfTargetPoints) const
2755 if(!arr || !arr->isAllocated())
2756 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationKriging::getValueOnMulti : input array is null or not allocated !");
2757 int nbOfRows(getNumberOfMeshPlaces(mesh));
2758 if(arr->getNumberOfTuples()!=nbOfRows)
2760 std::ostringstream oss; oss << "MEDCouplingFieldDiscretizationKriging::getValueOnMulti : input array does not have correct number of tuples ! Excepted " << nbOfRows << " having " << arr->getNumberOfTuples() << " !";
2761 throw INTERP_KERNEL::Exception(oss.str().c_str());
2763 int nbCols(-1),nbCompo(arr->getNumberOfComponents());
2764 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> m(computeEvaluationMatrixOnGivenPts(mesh,loc,nbOfTargetPoints,nbCols));
2765 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret(DataArrayDouble::New());
2766 ret->alloc(nbOfTargetPoints,nbCompo);
2767 INTERP_KERNEL::matrixProduct(m->begin(),nbOfTargetPoints,nbCols,arr->begin(),nbOfRows,nbCompo,ret->getPointer());
2771 void MEDCouplingFieldDiscretizationKriging::reprQuickOverview(std::ostream& stream) const
2773 stream << "Kriging spatial discretization.";
2777 * Returns the matrix of size nbRows = \a nbOfTargetPoints and \a nbCols = \a nbCols. This matrix is useful if
2779 * \return the new result matrix to be deallocated by the caller.
2781 DataArrayDouble *MEDCouplingFieldDiscretizationKriging::computeEvaluationMatrixOnGivenPts(const MEDCouplingMesh *mesh, const double *loc, int nbOfTargetPoints, int& nbCols) const
2783 int isDrift(-1),nbRows(-1);
2784 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrixInv(computeInverseMatrix(mesh,isDrift,nbRows));
2786 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coords=getLocalizationOfDiscValues(mesh);
2787 int nbOfPts(coords->getNumberOfTuples()),dimension(coords->getNumberOfComponents());
2788 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> locArr=DataArrayDouble::New();
2789 locArr->useArray(loc,false,CPP_DEALLOC,nbOfTargetPoints,dimension);
2792 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrix2=coords->buildEuclidianDistanceDenseMatrixWith(locArr);
2793 operateOnDenseMatrix(mesh->getSpaceDimension(),nbOfTargetPoints*nbOfPts,matrix2->getPointer());
2795 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrix3=DataArrayDouble::New();
2796 matrix3->alloc(nbOfTargetPoints*nbRows,1);
2797 double *work=matrix3->getPointer();
2798 const double *workCst(matrix2->begin()),*workCst2(loc);
2799 for(int i=0;i<nbOfTargetPoints;i++,workCst+=nbOfPts,workCst2+=isDrift-1)
2801 for(int j=0;j<nbOfPts;j++)
2802 work[i*nbRows+j]=workCst[j];
2803 work[i*nbRows+nbOfPts]=1.0;
2804 for(int j=0;j<isDrift-1;j++)
2805 work[i*nbRows+(nbOfPts+1+j)]=workCst2[j];
2807 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret(DataArrayDouble::New());
2808 ret->alloc(nbOfTargetPoints,nbRows);
2809 INTERP_KERNEL::matrixProduct(matrix3->begin(),nbOfTargetPoints,nbRows,matrixInv->begin(),nbRows,nbRows,ret->getPointer());
2810 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret2(DataArrayDouble::New());
2811 ret2->alloc(nbOfTargetPoints*nbOfPts,1);
2812 workCst=ret->begin(); work=ret2->getPointer();
2813 for(int i=0;i<nbOfTargetPoints;i++,workCst+=nbRows)
2814 work=std::copy(workCst,workCst+nbOfPts,work);
2819 * 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
2820 * when multiplied by the vector of values attached to each point.
2822 * \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.
2823 * \param [out] matSz the size of returned square matrix
2824 * \return the new result matrix to be deallocated by the caller.
2826 DataArrayDouble *MEDCouplingFieldDiscretizationKriging::computeInverseMatrix(const MEDCouplingMesh *mesh, int& isDrift, int& matSz) const
2829 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationKriging::computeVectorOfCoefficients : NULL input mesh !");
2830 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coords=getLocalizationOfDiscValues(mesh);
2831 int nbOfPts=coords->getNumberOfTuples();
2832 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrix=coords->buildEuclidianDistanceDenseMatrix();
2833 operateOnDenseMatrix(mesh->getSpaceDimension(),nbOfPts*nbOfPts,matrix->getPointer());
2835 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrixWithDrift=performDrift(matrix,coords,isDrift);
2836 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrixInv=DataArrayDouble::New();
2837 matSz=nbOfPts+isDrift;
2838 matrixInv->alloc(matSz*matSz,1);
2839 INTERP_KERNEL::inverseMatrix(matrixWithDrift->getConstPointer(),matSz,matrixInv->getPointer());
2840 return matrixInv.retn();
2844 * 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
2845 * number of tuples should be equal to the number of representing points in \a mesh.
2847 * \param [in] mesh is the sources of nodes on which kriging will be done regarding the parameters and the value of \c this->getSpaceDimension()
2848 * \param [in] arr input field DataArrayDouble whose number of tuples must be equal to the number of nodes in \a mesh
2849 * \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.
2850 * Whatever the value of \a isDrift the number of tuples of returned DataArrayDouble will be equal to \c arr->getNumberOfTuples() + \a isDrift.
2851 * \return a newly allocated array containing coefficients including or not drift coefficient at the end depending the value of \a isDrift parameter.
2853 DataArrayDouble *MEDCouplingFieldDiscretizationKriging::computeVectorOfCoefficients(const MEDCouplingMesh *mesh, const DataArrayDouble *arr, int& isDrift) const
2856 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> matrixInv(computeInverseMatrix(mesh,isDrift,nbRows));
2857 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> KnewiK=DataArrayDouble::New();
2858 KnewiK->alloc(nbRows*1,1);
2859 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arr2=DataArrayDouble::New();
2860 arr2->alloc(nbRows*1,1);
2861 double *work=std::copy(arr->begin(),arr->end(),arr2->getPointer());
2862 std::fill(work,work+isDrift,0.);
2863 INTERP_KERNEL::matrixProduct(matrixInv->getConstPointer(),nbRows,nbRows,arr2->getConstPointer(),nbRows,1,KnewiK->getPointer());
2864 return KnewiK.retn();
2868 * 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.
2870 * \param [in] spaceDimension space dimension of the input mesh on which the Kriging has to be performed
2871 * \param [in] nbOfElems is the result of the product of nb of rows and the nb of columns of matrix \a matrixPtr
2872 * \param [in,out] matrixPtr is the dense matrix whose on each values the operation will be applied
2874 void MEDCouplingFieldDiscretizationKriging::operateOnDenseMatrix(int spaceDimension, int nbOfElems, double *matrixPtr) const
2876 switch(spaceDimension)
2880 for(int i=0;i<nbOfElems;i++)
2882 double val=matrixPtr[i];
2883 matrixPtr[i]=val*val*val;
2889 for(int i=0;i<nbOfElems;i++)
2891 double val=matrixPtr[i];
2893 matrixPtr[i]=val*val*log(val);
2899 //nothing here : it is not a bug g(h)=h with spaceDim 3.
2903 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationKriging::operateOnDenseMatrix : only dimension 1, 2 and 3 implemented !");
2908 * Starting from a square matrix \a matr, this method returns a newly allocated dense square matrix whose \a matr is included in returned matrix
2909 * in the top left corner, and in the remaining returned matrix the parameters to take into account about the kriging drift.
2910 * For the moment only linear srift is implemented.
2912 * \param [in] arr the position of points were input mesh geometry is considered for Kriging
2913 * \param [in] matr input matrix whose drift part will be added
2914 * \param [out] delta the difference between the size of the output matrix and the input matrix \a matr.
2915 * \return a newly allocated matrix bigger than input matrix \a matr.
2917 DataArrayDouble *MEDCouplingFieldDiscretizationKriging::performDrift(const DataArrayDouble *matr, const DataArrayDouble *arr, int& delta) const
2919 int spaceDimension=arr->getNumberOfComponents();
2920 delta=spaceDimension+1;
2921 int szOfMatrix=arr->getNumberOfTuples();
2922 if(szOfMatrix*szOfMatrix!=matr->getNumberOfTuples())
2923 throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationKriging::performDrift : invalid size");
2924 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
2925 ret->alloc((szOfMatrix+delta)*(szOfMatrix+delta),1);
2926 const double *srcWork=matr->getConstPointer();
2927 const double *srcWork2=arr->getConstPointer();
2928 double *destWork=ret->getPointer();
2929 for(int i=0;i<szOfMatrix;i++)
2931 destWork=std::copy(srcWork,srcWork+szOfMatrix,destWork);
2932 srcWork+=szOfMatrix;
2934 destWork=std::copy(srcWork2,srcWork2+spaceDimension,destWork);
2935 srcWork2+=spaceDimension;
2937 std::fill(destWork,destWork+szOfMatrix,1.); destWork+=szOfMatrix;
2938 std::fill(destWork,destWork+spaceDimension+1,0.); destWork+=spaceDimension+1;
2939 MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arrNoI=arr->toNoInterlace();
2940 srcWork2=arrNoI->getConstPointer();
2941 for(int i=0;i<spaceDimension;i++)
2943 destWork=std::copy(srcWork2,srcWork2+szOfMatrix,destWork);
2944 srcWork2+=szOfMatrix;
2945 std::fill(destWork,destWork+spaceDimension+1,0.);
2946 destWork+=spaceDimension+1;