#include <shape_matrix.h>
Collaboration diagram for ShapeMatrix< _Vec_, Dim >:

Public Member Functions | |
| void | Allocate () |
| void | apply_changing (_Vec_ &v, double *multL, double *lambdas) |
| void | build_least_square_rhs (Vector *rhs, int at, int elem, double v) |
| void | build_node_shape (double *node_shape) |
| void | build_rhs (double *rhs, _Vec_ &v) |
| void | buildContinuConstraint (double *A, double *node_shape) |
| void | buildFormeFaibleCMatrix (FullMatrix &C, double *lambdas) |
| void | buildLeastSquareMatrix (FullMatrix *A) |
| void | changeAtomIndirection (int el, int old_ind_atome, int new_ind_atome) |
| void | correct_forme_faible (double &v, Vector &rhs, double *lamb, int at, unsigned int elem) |
| void | GetRelatedAtomsToElem (int elem, std::vector< unsigned int > &atoms) |
| void | GetRelatedNodesToElem (int elem, std::vector< unsigned int > &nodes) |
| int * | indirectionAtomTable () |
| double | InterpolateAtom (int at, _Vec_ &field, int elem, int deg) |
| void | operator() (int atome, int node_global, int node, int el, double value) |
| void | print () |
| void | SetSub (int i, int m, int n) |
| ShapeMatrix (int t) | |
| void | swapAtomIndirections () |
| ~ShapeMatrix () | |
Private Attributes | |
| int ** | currentIndirectionAtome |
| int * | indA_contigu [2] |
| int * | indGN_contigu |
| int ** | indirectionAtome [2] |
| int ** | indirectionGlobalNode |
| int ** | indirectionNode |
| int * | indN_contigu |
| int ** | newIndirectionAtome |
| int * | remplissageA |
| int * | remplissageN |
| FullMatrix ** | SubMatrix |
| double * | tab_contigu |
| int | taille |
| unsigned int * | tailles |
Definition at line 48 of file shape_matrix.h.
| ShapeMatrix< _Vec_, Dim >::ShapeMatrix | ( | int | t | ) | [inline] |
Definition at line 51 of file shape_matrix.h.
References ShapeMatrix< _Vec_, Dim >::currentIndirectionAtome, ShapeMatrix< _Vec_, Dim >::indirectionAtome, ShapeMatrix< _Vec_, Dim >::indirectionGlobalNode, ShapeMatrix< _Vec_, Dim >::indirectionNode, ShapeMatrix< _Vec_, Dim >::newIndirectionAtome, ShapeMatrix< _Vec_, Dim >::remplissageA, ShapeMatrix< _Vec_, Dim >::remplissageN, ShapeMatrix< _Vec_, Dim >::SubMatrix, ShapeMatrix< _Vec_, Dim >::tab_contigu, ShapeMatrix< _Vec_, Dim >::taille, and ShapeMatrix< _Vec_, Dim >::tailles.
00051 { 00052 SubMatrix = new FullMatrix *[t]; 00053 tailles = new unsigned int[2*t]; 00054 indirectionAtome[0] = new int *[t]; 00055 indirectionAtome[1] = new int *[t]; 00056 currentIndirectionAtome = indirectionAtome[0]; 00057 newIndirectionAtome = indirectionAtome[1]; 00058 indirectionNode = new int *[t]; 00059 indirectionGlobalNode = new int *[t]; 00060 remplissageA = new int[t]; 00061 remplissageN = new int[t]; 00062 memset(remplissageA,0,sizeof(int)*t); 00063 memset(remplissageN,0,sizeof(int)*t); 00064 00065 tab_contigu = NULL; 00066 taille = t; 00067 };
| ShapeMatrix< _Vec_, Dim >::~ShapeMatrix | ( | ) | [inline] |
Definition at line 69 of file shape_matrix.h.
References ShapeMatrix< _Vec_, Dim >::indA_contigu, ShapeMatrix< _Vec_, Dim >::indGN_contigu, ShapeMatrix< _Vec_, Dim >::indirectionAtome, ShapeMatrix< _Vec_, Dim >::indirectionGlobalNode, ShapeMatrix< _Vec_, Dim >::indirectionNode, ShapeMatrix< _Vec_, Dim >::indN_contigu, ShapeMatrix< _Vec_, Dim >::remplissageA, ShapeMatrix< _Vec_, Dim >::remplissageN, ShapeMatrix< _Vec_, Dim >::SubMatrix, ShapeMatrix< _Vec_, Dim >::tab_contigu, ShapeMatrix< _Vec_, Dim >::taille, and ShapeMatrix< _Vec_, Dim >::tailles.
00069 { 00070 for (int i=0;i<taille;++i) 00071 delete SubMatrix[i]; 00072 00073 if (tab_contigu) 00074 delete [] tab_contigu; 00075 if (indA_contigu[0]) 00076 delete [] indA_contigu[0]; 00077 if (indA_contigu[1]) 00078 delete [] indA_contigu[1]; 00079 if (indN_contigu) 00080 delete [] indN_contigu; 00081 if (indGN_contigu) 00082 delete [] indGN_contigu; 00083 00084 delete [] remplissageA; 00085 delete [] remplissageN; 00086 delete [] indirectionNode; 00087 delete [] indirectionGlobalNode; 00088 delete [] indirectionAtome[0]; 00089 delete [] indirectionAtome[1]; 00090 delete [] SubMatrix; 00091 delete [] tailles; 00092 }
| void ShapeMatrix< _Vec_, Dim >::Allocate | ( | ) | [inline] |
Definition at line 101 of file shape_matrix.h.
References ShapeMatrix< _Vec_, Dim >::indA_contigu, ShapeMatrix< _Vec_, Dim >::indGN_contigu, ShapeMatrix< _Vec_, Dim >::indirectionAtome, ShapeMatrix< _Vec_, Dim >::indirectionGlobalNode, ShapeMatrix< _Vec_, Dim >::indirectionNode, ShapeMatrix< _Vec_, Dim >::indN_contigu, ShapeMatrix< _Vec_, Dim >::SubMatrix, ShapeMatrix< _Vec_, Dim >::tab_contigu, ShapeMatrix< _Vec_, Dim >::taille, and ShapeMatrix< _Vec_, Dim >::tailles.
00101 { 00102 //premier passage pour calculer tout les 00103 //coeficients de la smatrice globale 00104 00105 int m,n; 00106 00107 int total_size=0; 00108 int total_atoms=0; 00109 int total_nodes=0; 00110 00111 for (int i=0;i<taille;++i){ 00112 m = tailles[i*2]; 00113 n = tailles[i*2+1]; 00114 00115 total_size+= (m+1)*n; 00116 total_atoms+=m; 00117 total_nodes+=n; 00118 } 00119 00120 //allocation 00121 tab_contigu = new double[total_size]; 00122 int index = 0; 00123 00124 indA_contigu[0] = new int[total_atoms]; 00125 indA_contigu[1] = new int[total_atoms]; 00126 int index_atoms = 0; 00127 00128 indN_contigu = new int[total_nodes]; 00129 indGN_contigu = new int[total_nodes]; 00130 int index_nodes = 0; 00131 00132 //deuxieme passage pour allouer les structures de control des matrices 00133 for (int i=0;i<taille;++i){ 00134 m = tailles[i*2]; 00135 n = tailles[i*2+1]; 00136 00137 SubMatrix[i] = new FullMatrix(m+1,n,tab_contigu+index); 00138 00139 indirectionAtome[0][i] = indA_contigu[0]+index_atoms; 00140 indirectionAtome[1][i] = indA_contigu[1]+index_atoms; 00141 indirectionNode[i] = indN_contigu+index_nodes; 00142 indirectionGlobalNode[i] = indGN_contigu+index_nodes; 00143 00144 //incrementation des indexs 00145 index += (m+1)*n; 00146 index_atoms += m; 00147 index_nodes += n; 00148 } 00149 }
| void ShapeMatrix< _Vec_, Dim >::apply_changing | ( | _Vec_ & | v, | |
| double * | multL, | |||
| double * | lambdas | |||
| ) | [inline] |
Definition at line 262 of file shape_matrix.h.
References ShapeMatrix< _Vec_, Dim >::currentIndirectionAtome, DBG_ALL, DUMP, ShapeMatrix< _Vec_, Dim >::indirectionGlobalNode, ShapeMatrix< _Vec_, Dim >::indirectionNode, ShapeMatrix< _Vec_, Dim >::SubMatrix, and ShapeMatrix< _Vec_, Dim >::taille.
00262 { 00263 // unsigned int nbT = v.size()/Dim; 00264 00265 // je parcours tout les blocs 00266 for (int k =0 ; k < taille ; ++k){ 00267 //je chope la sous matrice 00268 FullMatrix & mat = *SubMatrix[k]; 00269 //je chope les indirection sur les noeuds globaux comme 00270 // j'ai un acces direct au grand vecteur 00271 int * indGN = indirectionGlobalNode[k]; 00272 int * indN = indirectionNode[k]; 00273 int * indA = currentIndirectionAtome[k]; 00274 // je parcours tout les atomes du bloc 00275 for (int i = 0 ; i < mat.m()-1; ++i) 00276 //je parcours tout les noeuds 00277 for (int j = 0 ; j < mat.n(); ++j){ 00278 DUMP("adding v[" << indGN[j] << "] with " << -1.0*mat(i,j) << " * " << multL[indA[i]*Dim+1] << " {rhs[" << indA[i] << "]}/ " 00279 << lambdas[indN[j]]<< " = " << -1.0*mat(i,j)*multL[indA[i]*Dim+1]/lambdas[indN[j]],DBG_ALL); 00280 00281 v.add(indGN[j]*Dim,-1.0*mat(i,j)*multL[indA[i]*Dim]/lambdas[indN[j]]); 00282 if (Dim > 1) 00283 v.add(indGN[j]*Dim+1,-1.0*mat(i,j)*multL[indA[i]*Dim+1]/lambdas[indN[j]]); 00284 if (Dim == 3) 00285 v.add(indGN[j]*Dim+2,-1.0*mat(i,j)*multL[indA[i]*Dim+2]/lambdas[indN[j]]); 00286 } 00287 } 00288 };
| void ShapeMatrix< _Vec_, Dim >::build_least_square_rhs | ( | Vector * | rhs, | |
| int | at, | |||
| int | elem, | |||
| double | v | |||
| ) | [inline] |
Definition at line 379 of file shape_matrix.h.
References ShapeMatrix< _Vec_, Dim >::currentIndirectionAtome, ShapeMatrix< _Vec_, Dim >::indirectionNode, MyMatrix< T >::m(), MyMatrix< T >::n(), and ShapeMatrix< _Vec_, Dim >::SubMatrix.
00379 { 00380 int * indA = currentIndirectionAtome[elem]; 00381 int * indN = indirectionNode[elem]; 00382 FullMatrix & mat = *SubMatrix[elem]; 00383 00384 /* trouve l'atome concerné */ 00385 int i,j; 00386 for (i = 0 ; i < mat.m()-1; ++i) 00387 if (indA[i] == at) 00388 break; 00389 00390 /* construit le second membre de la methode least square */ 00391 00392 for (j = 0 ; j < mat.n(); ++j) 00393 (*rhs)(indN[j]) += mat(i,j)*v; 00394 00395 };
Here is the call graph for this function:

| void ShapeMatrix< _Vec_, Dim >::build_node_shape | ( | double * | node_shape | ) | [inline] |
Definition at line 202 of file shape_matrix.h.
References DBG_ALL, DUMP, ShapeMatrix< _Vec_, Dim >::indirectionNode, MyMatrix< T >::m(), MyMatrix< T >::n(), ShapeMatrix< _Vec_, Dim >::SubMatrix, and ShapeMatrix< _Vec_, Dim >::taille.
00202 { 00203 for (int k =0 ; k < taille ; ++k){ 00204 int * indN = indirectionNode[k]; 00205 FullMatrix & mat = *SubMatrix[k]; 00206 for (int j = 0 ; j < mat.n(); ++j){ 00207 DUMP("Summing node_shape[" << indN[j] << " (el= " << k << ", j = " << j << ")] with " << mat(mat.m()-1,j),DBG_ALL); 00208 node_shape[indN[j]] += mat(mat.m()-1,j); 00209 } 00210 } 00211 }
Here is the call graph for this function:

| void ShapeMatrix< _Vec_, Dim >::build_rhs | ( | double * | rhs, | |
| _Vec_ & | v | |||
| ) | [inline] |
Definition at line 227 of file shape_matrix.h.
References ShapeMatrix< _Vec_, Dim >::currentIndirectionAtome, DBG_ALL, DUMP, ShapeMatrix< _Vec_, Dim >::indirectionGlobalNode, ShapeMatrix< _Vec_, Dim >::SubMatrix, and ShapeMatrix< _Vec_, Dim >::taille.
00227 { 00228 00229 // unsigned int nbT = v.size()/Dim; 00230 00231 // je parcours tout les blocs 00232 for (int k =0 ; k < taille ; ++k){ 00233 //je chope la sous matrice 00234 FullMatrix & mat = *SubMatrix[k]; 00235 //je chope les indirection sur les noeuds globaux comme 00236 // j'ai un acces direct au grand vecteur 00237 int * indGN = indirectionGlobalNode[k]; 00238 int * indA = currentIndirectionAtome[k]; 00239 // je parcours tout les atomes du bloc 00240 for (int i = 0 ; i < mat.m()-1; ++i) 00241 //je parcours tout les noeuds 00242 for (int j = 0 ; j < mat.n(); ++j){ 00243 00244 DUMP("Before rhs[" << indA[i] << "]=" << rhs[Dim*indA[i]+1],DBG_ALL); 00245 00246 rhs[Dim*indA[i]] += v[indGN[j]*Dim] * mat(i,j); 00247 00248 if (Dim > 1) 00249 rhs[Dim*indA[i]+1] += v[indGN[j]*Dim+1] * mat(i,j); 00250 if (Dim ==3) 00251 rhs[Dim*indA[i]+2] += v[indGN[j]*Dim+2] * mat(i,j); 00252 00253 DUMP("build rhs to " << Dim*indA[i] << " " << Dim*indA[i]+ 1 << " " << Dim*indA[i]+3,DBG_ALL); 00254 DUMP("Building " << rhs[Dim*indA[i]+1] << "=rhs[" << indA[i] << "] with " << mat(i,j) << "*" << v[indGN[j]*Dim +1] << " = " 00255 << v[indGN[j]*Dim+1] * mat(i,j) << " (index = " << indGN[j] << ")",DBG_ALL); 00256 00257 } 00258 } 00259 }
| void ShapeMatrix< _Vec_, Dim >::buildContinuConstraint | ( | double * | A, | |
| double * | node_shape | |||
| ) | [inline] |
Definition at line 213 of file shape_matrix.h.
References ShapeMatrix< _Vec_, Dim >::currentIndirectionAtome, DBG_ALL, DUMP, ShapeMatrix< _Vec_, Dim >::indirectionNode, MyMatrix< T >::m(), MyMatrix< T >::n(), ShapeMatrix< _Vec_, Dim >::SubMatrix, and ShapeMatrix< _Vec_, Dim >::taille.
00213 { 00214 for (int k =0 ; k < taille ; ++k){ 00215 int * indA = currentIndirectionAtome[k]; 00216 int * indN = indirectionNode[k]; 00217 FullMatrix & mat = *SubMatrix[k]; 00218 for (int i = 0 ; i < mat.m()-1; ++i) 00219 for (int j = 0 ; j < mat.n(); ++j){ 00220 DUMP("Building A[" << indA[i] << "] with " << mat(i,j) << "*" << node_shape[indN[j]] << " = " 00221 << mat(i,j)*node_shape[indN[j]],DBG_ALL); 00222 A[indA[i]] += mat(i,j)*node_shape[indN[j]]; 00223 } 00224 } 00225 }
Here is the call graph for this function:

| void ShapeMatrix< _Vec_, Dim >::buildFormeFaibleCMatrix | ( | FullMatrix & | C, | |
| double * | lambdas | |||
| ) | [inline] |
Definition at line 400 of file shape_matrix.h.
References ShapeMatrix< _Vec_, Dim >::currentIndirectionAtome, DBG_ALL, DUMP, ShapeMatrix< _Vec_, Dim >::indirectionNode, MyMatrix< T >::m(), MyMatrix< T >::n(), ShapeMatrix< _Vec_, Dim >::SubMatrix, and ShapeMatrix< _Vec_, Dim >::taille.
00400 { 00401 00402 00403 /* le produit de deux shape sur la position d'un atome 00404 n'est nul que si les deux shapes font partie 00405 d'un element commun : on a une seul boucle pricipale sur 00406 les elements */ 00407 00408 00409 for (int k =0 ; k < taille ; ++k){ 00410 int * indN = indirectionNode[k]; 00411 int * indA = currentIndirectionAtome[k]; 00412 00413 FullMatrix & mat = *SubMatrix[k]; 00414 00415 /* ici on met la double boucle sur les noeuds pour calculer la sous matrice de la grande */ 00416 00417 for (int j1 = 0 ; j1 < mat.n() ; ++j1) 00418 for (int j2 = 0 ; j2 < mat.n() ; ++j2) 00419 { 00420 /* ici on met la boucle de la somme sur les atomes */ 00421 for (int i = 0 ; i < mat.m()-1 ; ++i){ 00422 00423 C(indN[j1],indN[j2]) += mat(i,j1) * mat(i,j2) / lambdas[indA[i]]; 00424 DUMP("C[" << indN[j1] << "," << indN[j2] << "] = " << C(indN[j1],indN[j2]),DBG_ALL); 00425 } 00426 } 00427 } 00428 00429 };
Here is the call graph for this function:

| void ShapeMatrix< _Vec_, Dim >::buildLeastSquareMatrix | ( | FullMatrix * | A | ) | [inline] |
Definition at line 345 of file shape_matrix.h.
References DBG_ALL, DUMP, ShapeMatrix< _Vec_, Dim >::indirectionNode, MyMatrix< T >::m(), MyMatrix< T >::n(), MyMatrix< T >::PrintToFile(), ShapeMatrix< _Vec_, Dim >::SubMatrix, ShapeMatrix< _Vec_, Dim >::taille, and FullMatrix::zero().
00345 { 00346 A->zero(); 00347 00348 /* on veut construire les A(h,k) */ 00349 00350 /* le produit de deux shape sur la position d'un atome 00351 n'est nul que si les deux shapes font partie 00352 d'un element commun : on a une seul boucle pricipale sur 00353 les elements */ 00354 00355 00356 for (int k =0 ; k < taille ; ++k){ 00357 int * indN = indirectionNode[k]; 00358 FullMatrix & mat = *SubMatrix[k]; 00359 00360 /* ici on met la double boucle sur les noeuds pour calculer la sous matrice de la grande */ 00361 00362 00363 for (int j1 = 0 ; j1 < mat.n() ; ++j1) 00364 for (int j2 = 0 ; j2 < mat.n() ; ++j2) 00365 { 00366 /* ici on met la boucle de la somme sur les atomes */ 00367 for (int i = 0 ; i < mat.m()-1 ; ++i){ 00368 (*A)(indN[j1],indN[j2]) += mat(i,j1) * mat(i,j2); 00369 DUMP("A[" << indN[j1] << "," << indN[j2] << "] = " << (*A)(indN[j1],indN[j2]),DBG_ALL); 00370 } 00371 } 00372 } 00373 00374 FILE * f = fopen("least-square-constraint.mat","wb+"); 00375 A->PrintToFile(f); 00376 fclose(f); 00377 }
Here is the call graph for this function:

| void ShapeMatrix< _Vec_, Dim >::changeAtomIndirection | ( | int | el, | |
| int | old_ind_atome, | |||
| int | new_ind_atome | |||
| ) | [inline] |
Definition at line 167 of file shape_matrix.h.
References ShapeMatrix< _Vec_, Dim >::currentIndirectionAtome, MyMatrix< T >::m(), ShapeMatrix< _Vec_, Dim >::newIndirectionAtome, and ShapeMatrix< _Vec_, Dim >::SubMatrix.
Referenced by BridgingPar< DomainA, DomainC, Dim, Pond >::UnfragmentIndirectionArrays().
00167 { 00168 FullMatrix & mat = *SubMatrix[el]; 00169 int * indA = currentIndirectionAtome[el]; 00170 int * indAN = newIndirectionAtome[el]; 00171 for (int i = 0 ; i < mat.m()-1;++i) 00172 if (indA[i] == old_ind_atome) { 00173 indAN[i] = new_ind_atome; 00174 break; 00175 } 00176 }
Here is the call graph for this function:

| void ShapeMatrix< _Vec_, Dim >::correct_forme_faible | ( | double & | v, | |
| Vector & | rhs, | |||
| double * | lamb, | |||
| int | at, | |||
| unsigned int | elem | |||
| ) | [inline] |
Definition at line 433 of file shape_matrix.h.
References ShapeMatrix< _Vec_, Dim >::currentIndirectionAtome, ShapeMatrix< _Vec_, Dim >::indirectionNode, MyMatrix< T >::m(), MyMatrix< T >::n(), and ShapeMatrix< _Vec_, Dim >::SubMatrix.
00433 { 00434 int * indN = indirectionNode[elem]; 00435 int * indA = currentIndirectionAtome[elem]; 00436 FullMatrix & mat = *SubMatrix[elem]; 00437 00438 /* trouve l'atome concerné */ 00439 int i,j; 00440 for (i = 0 ; i < mat.m()-1; ++i) 00441 if (indA[i] == at) 00442 break; 00443 00444 00445 for (j = 0 ; j < mat.n() ; ++j){ 00446 v -= mat(i,j)*rhs(indN[j])/lamb[at]; 00447 } 00448 00449 };
Here is the call graph for this function:

| void ShapeMatrix< _Vec_, Dim >::GetRelatedAtomsToElem | ( | int | elem, | |
| std::vector< unsigned int > & | atoms | |||
| ) | [inline] |
Definition at line 316 of file shape_matrix.h.
References ShapeMatrix< _Vec_, Dim >::currentIndirectionAtome, MyMatrix< T >::m(), and ShapeMatrix< _Vec_, Dim >::SubMatrix.
00316 { 00317 int * indA = currentIndirectionAtome[elem]; 00318 FullMatrix & mat = *SubMatrix[elem]; 00319 00320 atoms.clear(); 00321 00322 /* parcours les atomes */ 00323 int i; 00324 for (i = 0 ; i < mat.m()-1; ++i) 00325 atoms.push_back(indA[i]); 00326 };
Here is the call graph for this function:

| void ShapeMatrix< _Vec_, Dim >::GetRelatedNodesToElem | ( | int | elem, | |
| std::vector< unsigned int > & | nodes | |||
| ) | [inline] |
Definition at line 329 of file shape_matrix.h.
References ShapeMatrix< _Vec_, Dim >::indirectionGlobalNode, MyMatrix< T >::n(), and ShapeMatrix< _Vec_, Dim >::SubMatrix.
00329 { 00330 int * indGN = indirectionGlobalNode[elem]; 00331 FullMatrix & mat = *SubMatrix[elem]; 00332 00333 nodes.clear(); 00334 00335 /* parcours les noeuds */ 00336 int j; 00337 for (j = 0 ; j < mat.n(); ++j) 00338 nodes.push_back(indGN[j]); 00339 };
Here is the call graph for this function:

| int* ShapeMatrix< _Vec_, Dim >::indirectionAtomTable | ( | ) | [inline] |
Definition at line 554 of file shape_matrix.h.
References ShapeMatrix< _Vec_, Dim >::indA_contigu.
00554 { 00555 return indA_contigu; 00556 }
| double ShapeMatrix< _Vec_, Dim >::InterpolateAtom | ( | int | at, | |
| _Vec_ & | field, | |||
| int | elem, | |||
| int | deg | |||
| ) | [inline] |
Definition at line 290 of file shape_matrix.h.
References ShapeMatrix< _Vec_, Dim >::currentIndirectionAtome, DBG_ALL, DUMP, ShapeMatrix< _Vec_, Dim >::indirectionGlobalNode, MyMatrix< T >::m(), MyMatrix< T >::n(), and ShapeMatrix< _Vec_, Dim >::SubMatrix.
00290 { 00291 double res = 0; 00292 int * indA = currentIndirectionAtome[elem]; 00293 int * indGN = indirectionGlobalNode[elem]; 00294 FullMatrix & mat = *SubMatrix[elem]; 00295 00296 DUMP("interpolating for elem " << elem << " atom " << at,DBG_ALL); 00297 /* trouve l'atome concerné */ 00298 int i; 00299 for (i = 0 ; i < mat.m()-1; ++i) 00300 if (indA[i] == at) 00301 break; 00302 00303 /* effectue l'interpolation */ 00304 00305 // unsigned int nbT = field.size()/Dim; 00306 00307 for (int j = 0 ; j < mat.n(); ++j){ 00308 res += mat(i,j)*field[indGN[j]*Dim+deg]; 00309 } 00310 00311 return res; 00312 00313 };
Here is the call graph for this function:

| void ShapeMatrix< _Vec_, Dim >::operator() | ( | int | atome, | |
| int | node_global, | |||
| int | node, | |||
| int | el, | |||
| double | value | |||
| ) | [inline] |
Definition at line 178 of file shape_matrix.h.
References ShapeMatrix< _Vec_, Dim >::currentIndirectionAtome, DBG_ALL, DUMP, ShapeMatrix< _Vec_, Dim >::indirectionGlobalNode, ShapeMatrix< _Vec_, Dim >::indirectionNode, MyMatrix< T >::m(), MyMatrix< T >::n(), ShapeMatrix< _Vec_, Dim >::remplissageA, ShapeMatrix< _Vec_, Dim >::remplissageN, and ShapeMatrix< _Vec_, Dim >::SubMatrix.
00178 { 00179 int * indA = currentIndirectionAtome[el]; 00180 indA[remplissageA[el]]=atome; 00181 00182 int * indN = indirectionNode[el]; 00183 indN[remplissageN[el]]=node; 00184 00185 DUMP("setting shapes for atome " << atome << " and element " << el << " node " << node << " node_global " << node_global << " value " << value,DBG_ALL); 00186 int * indGN = indirectionGlobalNode[el]; 00187 indGN[remplissageN[el]]=node_global; 00188 00189 // je prend la matrice dense qui correspond a mon element 00190 FullMatrix & m = *SubMatrix[el]; 00191 m(remplissageA[el],remplissageN[el]) = value; 00192 // fait les sommes partielles sur les noeuds 00193 m(m.m()-1,remplissageN[el])+= value; 00194 00195 ++remplissageN[el]; 00196 if (remplissageN[el] == m.n()){ 00197 remplissageN[el]=0; 00198 ++remplissageA[el]; 00199 } 00200 }
Here is the call graph for this function:

| void ShapeMatrix< _Vec_, Dim >::print | ( | ) | [inline] |
Definition at line 151 of file shape_matrix.h.
References DBG_DETAIL, DUMP, ShapeMatrix< _Vec_, Dim >::indirectionAtome, MyMatrix< T >::m(), MyMatrix< T >::print(), ShapeMatrix< _Vec_, Dim >::remplissageA, ShapeMatrix< _Vec_, Dim >::SubMatrix, and ShapeMatrix< _Vec_, Dim >::taille.
00151 { 00152 for (int k=0;k<taille;++k){ 00153 if (!remplissageA[k]) return; 00154 FullMatrix & mat = *SubMatrix[k]; 00155 00156 DUMP("affiche les infos de " << k,DBG_DETAIL); 00157 00158 for (int i =0; i < mat.m() -1; ++i) 00159 std::cout << indirectionAtome[k][i] << "\t"; 00160 00161 std::cout << endl; 00162 00163 SubMatrix[k]->print(); 00164 } 00165 }
Here is the call graph for this function:

| void ShapeMatrix< _Vec_, Dim >::SetSub | ( | int | i, | |
| int | m, | |||
| int | n | |||
| ) | [inline] |
Definition at line 94 of file shape_matrix.h.
References FATAL, and ShapeMatrix< _Vec_, Dim >::tailles.
00094 { 00095 if (m <= 0 || n <= 0) 00096 FATAL("can not set a smatrix bloc size of zero : m=" << m << " n=" <<n); 00097 tailles[i*2] = m; 00098 tailles[i*2+1] = n; 00099 };
| void ShapeMatrix< _Vec_, Dim >::swapAtomIndirections | ( | ) | [inline] |
Definition at line 558 of file shape_matrix.h.
References ShapeMatrix< _Vec_, Dim >::currentIndirectionAtome, and ShapeMatrix< _Vec_, Dim >::newIndirectionAtome.
Referenced by BridgingPar< DomainA, DomainC, Dim, Pond >::UnfragmentIndirectionArrays().
00558 { 00559 int ** temp = currentIndirectionAtome; 00560 currentIndirectionAtome = newIndirectionAtome; 00561 newIndirectionAtome = temp; 00562 }
int** ShapeMatrix< _Vec_, Dim >::currentIndirectionAtome [private] |
Definition at line 573 of file shape_matrix.h.
Referenced by ShapeMatrix< _Vec_, Dim >::apply_changing(), ShapeMatrix< _Vec_, Dim >::build_least_square_rhs(), ShapeMatrix< _Vec_, Dim >::build_rhs(), ShapeMatrix< _Vec_, Dim >::buildContinuConstraint(), ShapeMatrix< _Vec_, Dim >::buildFormeFaibleCMatrix(), ShapeMatrix< _Vec_, Dim >::changeAtomIndirection(), ShapeMatrix< _Vec_, Dim >::correct_forme_faible(), ShapeMatrix< _Vec_, Dim >::GetRelatedAtomsToElem(), ShapeMatrix< _Vec_, Dim >::InterpolateAtom(), ShapeMatrix< _Vec_, Dim >::operator()(), ShapeMatrix< _Vec_, Dim >::ShapeMatrix(), and ShapeMatrix< _Vec_, Dim >::swapAtomIndirections().
int* ShapeMatrix< _Vec_, Dim >::indA_contigu[2] [private] |
Definition at line 587 of file shape_matrix.h.
Referenced by ShapeMatrix< _Vec_, Dim >::Allocate(), ShapeMatrix< _Vec_, Dim >::indirectionAtomTable(), and ShapeMatrix< _Vec_, Dim >::~ShapeMatrix().
int* ShapeMatrix< _Vec_, Dim >::indGN_contigu [private] |
Definition at line 590 of file shape_matrix.h.
Referenced by ShapeMatrix< _Vec_, Dim >::Allocate(), and ShapeMatrix< _Vec_, Dim >::~ShapeMatrix().
int** ShapeMatrix< _Vec_, Dim >::indirectionAtome[2] [private] |
Definition at line 576 of file shape_matrix.h.
Referenced by ShapeMatrix< _Vec_, Dim >::Allocate(), ShapeMatrix< _Vec_, Dim >::print(), ShapeMatrix< _Vec_, Dim >::ShapeMatrix(), and ShapeMatrix< _Vec_, Dim >::~ShapeMatrix().
int** ShapeMatrix< _Vec_, Dim >::indirectionGlobalNode [private] |
Definition at line 579 of file shape_matrix.h.
Referenced by ShapeMatrix< _Vec_, Dim >::Allocate(), ShapeMatrix< _Vec_, Dim >::apply_changing(), ShapeMatrix< _Vec_, Dim >::build_rhs(), ShapeMatrix< _Vec_, Dim >::GetRelatedNodesToElem(), ShapeMatrix< _Vec_, Dim >::InterpolateAtom(), ShapeMatrix< _Vec_, Dim >::operator()(), ShapeMatrix< _Vec_, Dim >::ShapeMatrix(), and ShapeMatrix< _Vec_, Dim >::~ShapeMatrix().
int** ShapeMatrix< _Vec_, Dim >::indirectionNode [private] |
Definition at line 578 of file shape_matrix.h.
Referenced by ShapeMatrix< _Vec_, Dim >::Allocate(), ShapeMatrix< _Vec_, Dim >::apply_changing(), ShapeMatrix< _Vec_, Dim >::build_least_square_rhs(), ShapeMatrix< _Vec_, Dim >::build_node_shape(), ShapeMatrix< _Vec_, Dim >::buildContinuConstraint(), ShapeMatrix< _Vec_, Dim >::buildFormeFaibleCMatrix(), ShapeMatrix< _Vec_, Dim >::buildLeastSquareMatrix(), ShapeMatrix< _Vec_, Dim >::correct_forme_faible(), ShapeMatrix< _Vec_, Dim >::operator()(), ShapeMatrix< _Vec_, Dim >::ShapeMatrix(), and ShapeMatrix< _Vec_, Dim >::~ShapeMatrix().
int* ShapeMatrix< _Vec_, Dim >::indN_contigu [private] |
Definition at line 589 of file shape_matrix.h.
Referenced by ShapeMatrix< _Vec_, Dim >::Allocate(), and ShapeMatrix< _Vec_, Dim >::~ShapeMatrix().
int** ShapeMatrix< _Vec_, Dim >::newIndirectionAtome [private] |
Definition at line 574 of file shape_matrix.h.
Referenced by ShapeMatrix< _Vec_, Dim >::changeAtomIndirection(), ShapeMatrix< _Vec_, Dim >::ShapeMatrix(), and ShapeMatrix< _Vec_, Dim >::swapAtomIndirections().
int* ShapeMatrix< _Vec_, Dim >::remplissageA [private] |
Definition at line 580 of file shape_matrix.h.
Referenced by ShapeMatrix< _Vec_, Dim >::operator()(), ShapeMatrix< _Vec_, Dim >::print(), ShapeMatrix< _Vec_, Dim >::ShapeMatrix(), and ShapeMatrix< _Vec_, Dim >::~ShapeMatrix().
int* ShapeMatrix< _Vec_, Dim >::remplissageN [private] |
Definition at line 581 of file shape_matrix.h.
Referenced by ShapeMatrix< _Vec_, Dim >::operator()(), ShapeMatrix< _Vec_, Dim >::ShapeMatrix(), and ShapeMatrix< _Vec_, Dim >::~ShapeMatrix().
FullMatrix** ShapeMatrix< _Vec_, Dim >::SubMatrix [private] |
Definition at line 572 of file shape_matrix.h.
Referenced by ShapeMatrix< _Vec_, Dim >::Allocate(), ShapeMatrix< _Vec_, Dim >::apply_changing(), ShapeMatrix< _Vec_, Dim >::build_least_square_rhs(), ShapeMatrix< _Vec_, Dim >::build_node_shape(), ShapeMatrix< _Vec_, Dim >::build_rhs(), ShapeMatrix< _Vec_, Dim >::buildContinuConstraint(), ShapeMatrix< _Vec_, Dim >::buildFormeFaibleCMatrix(), ShapeMatrix< _Vec_, Dim >::buildLeastSquareMatrix(), ShapeMatrix< _Vec_, Dim >::changeAtomIndirection(), ShapeMatrix< _Vec_, Dim >::correct_forme_faible(), ShapeMatrix< _Vec_, Dim >::GetRelatedAtomsToElem(), ShapeMatrix< _Vec_, Dim >::GetRelatedNodesToElem(), ShapeMatrix< _Vec_, Dim >::InterpolateAtom(), ShapeMatrix< _Vec_, Dim >::operator()(), ShapeMatrix< _Vec_, Dim >::print(), ShapeMatrix< _Vec_, Dim >::ShapeMatrix(), and ShapeMatrix< _Vec_, Dim >::~ShapeMatrix().
double* ShapeMatrix< _Vec_, Dim >::tab_contigu [private] |
Definition at line 586 of file shape_matrix.h.
Referenced by ShapeMatrix< _Vec_, Dim >::Allocate(), ShapeMatrix< _Vec_, Dim >::ShapeMatrix(), and ShapeMatrix< _Vec_, Dim >::~ShapeMatrix().
int ShapeMatrix< _Vec_, Dim >::taille [private] |
Definition at line 568 of file shape_matrix.h.
Referenced by ShapeMatrix< _Vec_, Dim >::Allocate(), ShapeMatrix< _Vec_, Dim >::apply_changing(), ShapeMatrix< _Vec_, Dim >::build_node_shape(), ShapeMatrix< _Vec_, Dim >::build_rhs(), ShapeMatrix< _Vec_, Dim >::buildContinuConstraint(), ShapeMatrix< _Vec_, Dim >::buildFormeFaibleCMatrix(), ShapeMatrix< _Vec_, Dim >::buildLeastSquareMatrix(), ShapeMatrix< _Vec_, Dim >::print(), ShapeMatrix< _Vec_, Dim >::ShapeMatrix(), and ShapeMatrix< _Vec_, Dim >::~ShapeMatrix().
unsigned int* ShapeMatrix< _Vec_, Dim >::tailles [private] |
Definition at line 570 of file shape_matrix.h.
Referenced by ShapeMatrix< _Vec_, Dim >::Allocate(), ShapeMatrix< _Vec_, Dim >::SetSub(), ShapeMatrix< _Vec_, Dim >::ShapeMatrix(), and ShapeMatrix< _Vec_, Dim >::~ShapeMatrix().
1.5.2