Changeset 4229
- Timestamp:
- 06/25/10 09:23:11 (15 years ago)
- Location:
- issm/trunk/src/c
- Files:
-
- 6 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk/src/c/DataSet/DataSet.h
r4221 r4229 130 130 /*}}}*/ 131 131 /*numerics: {{{1*/ 132 void DistributeDofs(int numberofnodes,int numdofspernode );133 void FlagClones(int numberofnodes );132 void DistributeDofs(int numberofnodes,int numdofspernode,int analysis_type); 133 void FlagClones(int numberofnodes,int analysis_type); 134 134 void FlagNodeSets(Vec pv_g, Vec pv_m, Vec pv_n, Vec pv_f, Vec pv_s,int analysis_type); 135 135 int NumberOfDofs(int analysis_type); 136 int NumberOfNodes(int analysis_type); 136 137 int NumberOfNodes(void); 137 void Ranks(int* ranks );138 void Ranks(int* ranks,int analysis_type); 138 139 /*}}}*/ 139 140 -
issm/trunk/src/c/DataSet/Nodes.cpp
r4221 r4229 44 44 /*Numerics*/ 45 45 /*FUNCTION Nodes::DistributeDofs{{{1*/ 46 void Nodes::DistributeDofs(int numberofobjects,int numberofdofsperobject ){46 void Nodes::DistributeDofs(int numberofobjects,int numberofdofsperobject,int analysis_type){ 47 47 48 48 extern int num_procs; … … 59 59 for (i=0;i<this->Size();i++){ 60 60 Node* node=(Node*)this->GetObjectByOffset(i); 61 node->DistributeDofs(&dofcount); 61 62 /*Check that this node corresponds to our analysis currently being carried out: */ 63 if (node->InAnalysis(analysis_type)){ 64 node->DistributeDofs(&dofcount); 65 } 62 66 } 63 67 … … 85 89 /*Ok, now every cpu knows where his dofs should start. Update the dof count: */ 86 90 for (i=0;i<this->Size();i++){ 87 Node* node=(Node*)this->GetObjectByOffset(i); 88 node->OffsetDofs(dofcount); 91 92 /*Check that this node corresponds to our analysis currently being carried out: */ 93 Node* node=(Node*)this->GetObjectByOffset(i); 94 if (node->InAnalysis(analysis_type)){ 95 node->OffsetDofs(dofcount); 96 } 97 89 98 } 90 99 … … 96 105 97 106 for (i=0;i<this->Size();i++){ 98 Node* node=(Node*)this->GetObjectByOffset(i); 99 node->ShowTrueDofs(truedofs); 107 /*Check that this node corresponds to our analysis currently being carried out: */ 108 Node* node=(Node*)this->GetObjectByOffset(i); 109 if (node->InAnalysis(analysis_type)){ 110 node->ShowTrueDofs(truedofs); 111 } 100 112 } 101 113 … … 104 116 /*Ok, now every cpu knows the true dofs of everyone else that is not a clone. Let the clones recover those true dofs: */ 105 117 for (i=0;i<this->Size();i++){ 106 Node* node=(Node*)this->GetObjectByOffset(i); 107 node->UpdateCloneDofs(alltruedofs); 118 /*Check that this node corresponds to our analysis currently being carried out: */ 119 Node* node=(Node*)this->GetObjectByOffset(i); 120 if (node->InAnalysis(analysis_type)){ 121 node->UpdateCloneDofs(alltruedofs); 122 } 108 123 } 109 124 … … 116 131 /*}}}*/ 117 132 /*FUNCTION Nodes::FlagClones{{{1*/ 118 void Nodes::FlagClones(int numberofobjects ){133 void Nodes::FlagClones(int numberofobjects,int analysis_type){ 119 134 120 135 int i; … … 131 146 132 147 /*Now go through all our objects and ask them to report to who they belong (which rank): */ 133 Ranks(ranks );148 Ranks(ranks,analysis_type); 134 149 135 150 /*We need to take the minimum rank for each vertex, and every cpu needs to get that result. That way, … … 141 156 /*Now go through all objects, and use minranks to flag which objects are cloned: */ 142 157 for(i=0;i<this->Size();i++){ 143 /*For this object, decide whether it is a clone: */ 144 Node* node=(Node*)this->GetObjectByOffset(i); 145 node->SetClone(minranks); 158 159 Node* node=(Node*)this->GetObjectByOffset(i); 160 161 /*Check that this node corresponds to our analysis currently being carried out: */ 162 if (node->InAnalysis(analysis_type)){ 163 164 /*For this object, decide whether it is a clone: */ 165 node->SetClone(minranks); 166 } 146 167 } 147 168 … … 158 179 159 180 for(i=0;i<this->Size();i++){ 181 160 182 Node* node=(Node*)this->GetObjectByOffset(i); 161 183 162 if (node->InAnalysis(analysis_type)){ 184 /*Check that this node corresponds to our analysis currently being carried out: */ 185 if (node->InAnalysis(analysis_type)){ 186 163 187 /*Plug set values intp our 4 set vectors: */ 164 188 node->CreateVecSets(pv_g,pv_m,pv_n,pv_f,pv_s); … … 215 239 } 216 240 /*}}}*/ 217 /*FUNCTION Nodes::NumberOfNodes {{{1*/241 /*FUNCTION Nodes::NumberOfNodes(){{{1*/ 218 242 int Nodes::NumberOfNodes(void){ 219 243 … … 244 268 } 245 269 /*}}}*/ 270 /*FUNCTION Nodes::NumberOfNodes{{{1*/ 271 int Nodes::NumberOfNodes(int analysis_type){ 272 273 int i; 274 275 int max_sid=0; 276 int sid; 277 int node_max_sid; 278 279 for(i=0;i<this->Size();i++){ 280 281 Node* node=(Node*)this->GetObjectByOffset(i); 282 283 /*Check that this node corresponds to our analysis currently being carried out: */ 284 if (node->InAnalysis(analysis_type)){ 285 286 sid=node->Sid(); 287 if (sid>max_sid)max_sid=sid; 288 } 289 } 290 291 #ifdef _PARALLEL_ 292 MPI_Reduce (&max_sid,&node_max_sid,1,MPI_INT,MPI_MAX,0,MPI_COMM_WORLD ); 293 MPI_Bcast(&node_max_sid,1,MPI_INT,0,MPI_COMM_WORLD); 294 max_sid=node_max_sid; 295 #endif 296 297 /*sid starts at 0*/ 298 max_sid++; 299 300 /*return*/ 301 return max_sid; 302 } 303 /*}}}*/ 246 304 /*FUNCTION Nodes::Ranks{{{1*/ 247 void Nodes::Ranks(int* ranks ){305 void Nodes::Ranks(int* ranks,int analysis_type){ 248 306 249 307 /*Go through nodes, and for each object, report it cpu: */ … … 254 312 255 313 for(i=0;i<this->Size();i++){ 256 Node* node=(Node*)this->GetObjectByOffset(i); 257 rank=node->MyRank(); 258 sid=node->Sid(); 259 260 /*Plug rank into ranks, according to id: */ 261 ranks[sid]=rank; 314 315 Node* node=(Node*)this->GetObjectByOffset(i); 316 317 /*Check that this node corresponds to our analysis currently being carried out: */ 318 if (node->InAnalysis(analysis_type)){ 319 320 rank=node->MyRank(); 321 sid=node->Sid(); 322 323 /*Plug rank into ranks, according to id: */ 324 ranks[sid]=rank; 325 } 262 326 } 263 327 return; -
issm/trunk/src/c/modules/BuildNodeSetsx/BuildNodeSetsx.cpp
r4211 r4229 77 77 VecFree(&vec_pv_s); 78 78 79 80 79 /*Create NodeSets* object: */ 81 80 nodesets=new NodeSets(pv_m,pv_n,pv_f,pv_s,gsize,msize,nsize,fsize,ssize); -
issm/trunk/src/c/modules/NodesDofx/NodesDofx.cpp
r4211 r4229 10 10 #include "../../EnumDefinitions/EnumDefinitions.h" 11 11 12 void NodesDofx(Nodes* nodes, Parameters* parameters ){12 void NodesDofx(Nodes* nodes, Parameters* parameters,int analysis_type){ 13 13 14 14 int noerr=1; … … 20 20 int numberofdofspernode; 21 21 22 /*First, recover number of vertices and nodes: */23 numberofnodes=nodes->NumberOfNodes( );22 /*First, recover number nodes associated th this analysis: */ 23 numberofnodes=nodes->NumberOfNodes(analysis_type); 24 24 25 25 /*Recover number of dofs per node: */ 26 26 found=parameters->FindParam(&numberofdofspernode,NumberOfDofsPerNodeEnum); 27 if(!found)ISSMERROR(" 27 if(!found)ISSMERROR("could not find numberofdofspernode in parameters"); 28 28 29 29 /*Ensure that only for each cpu, the partition border nodes only will be taken into account once 30 30 * across the cluster. To do so, we flag all the clone nodes: */ 31 nodes->FlagClones(numberofnodes );31 nodes->FlagClones(numberofnodes,analysis_type); 32 32 33 33 /*Go through all nodes, and build degree of freedom lists. Each node gets a fixed number of dofs. When 34 34 *a node has already been distributed dofs on one cpu, all other cpus with the same node cannot distribute it 35 35 *anymore. Use clone field to be sure of that: */ 36 nodes->DistributeDofs(numberofnodes,numberofdofspernode );36 nodes->DistributeDofs(numberofnodes,numberofdofspernode,analysis_type); 37 37 38 38 } -
issm/trunk/src/c/modules/NodesDofx/NodesDofx.h
r4211 r4229 10 10 11 11 /* local prototypes: */ 12 void NodesDofx(Nodes* nodes, Parameters* parameters );12 void NodesDofx(Nodes* nodes, Parameters* parameters,int analysis_type); 13 13 14 14 #endif /* _NODESDOFX_H */ -
issm/trunk/src/c/objects/FemModel.cpp
r4181 r4229 63 63 _printf_(" create degrees of freedom\n"); 64 64 VerticesDofx(&partition,&tpartition,vertices,parameters); 65 NodesDofx(nodes,parameters );65 NodesDofx(nodes,parameters,analysis_type); 66 66 67 67 _printf_(" create single point constraints\n");
Note:
See TracChangeset
for help on using the changeset viewer.