/[escript]/trunk/finley/src/Mesh_createNodeFileMappings.c
ViewVC logotype

Contents of /trunk/finley/src/Mesh_createNodeFileMappings.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1770 - (show annotations)
Mon Sep 8 06:58:47 2008 UTC (11 years ago) by ksteube
File MIME type: text/plain
File size: 17934 byte(s)
Detects an array reference out of bounds when running run_inputOutput.py

1
2 /* $Id:$ */
3
4 /*******************************************************
5 *
6 * Copyright 2003-2007 by ACceSS MNRF
7 * Copyright 2007 by University of Queensland
8 *
9 * http://esscc.uq.edu.au
10 * Primary Business: Queensland, Australia
11 * Licensed under the Open Software License version 3.0
12 * http://www.opensource.org/licenses/osl-3.0.php
13 *
14 *******************************************************/
15
16 /**************************************************************/
17
18 /* Finley: NodeFile : creates the mappings using the indexReducedNodes */
19 /* no distribution is happening */
20
21 /**************************************************************/
22
23 #include "Mesh.h"
24 #define UNUSED -1
25
26 /**************************************************************/
27
28 void Mesh_createDOFMappingAndCoupling(Finley_Mesh* in, bool_t use_reduced_elements)
29 {
30 index_t min_DOF, max_DOF, *shared=NULL, *offsetInShared=NULL, *locDOFMask=NULL, i, k, myFirstDOF, myLastDOF, *nodeMask=NULL, firstDOF, lastDOF, *globalDOFIndex, *wanted_DOFs=NULL;
31 dim_t mpiSize, len_loc_dof, numNeighbors, n, lastn, numNodes,*rcv_len=NULL, *snd_len=NULL, count;
32 Paso_MPI_rank myRank,p,p_min,p_max, *neighbor=NULL;
33 Paso_SharedComponents *rcv_shcomp=NULL, *snd_shcomp=NULL;
34 Finley_NodeMapping *this_mapping=NULL;
35 Paso_Connector* this_connector=NULL;
36 Paso_Distribution* dof_distribution;
37 Paso_MPIInfo *mpi_info = in->MPIInfo;
38 #ifdef PASO_MPI
39 MPI_Request* mpi_requests=NULL;
40 MPI_Status* mpi_stati=NULL;
41 #else
42 int *mpi_requests=NULL, *mpi_stati=NULL;
43 #endif
44
45 numNodes=in->Nodes->numNodes;
46 if (use_reduced_elements) {
47 dof_distribution=in->Nodes->reducedDegreesOfFreedomDistribution;
48 globalDOFIndex=in->Nodes->globalReducedDOFIndex;
49 } else {
50 dof_distribution=in->Nodes->degreesOfFreedomDistribution;
51 globalDOFIndex=in->Nodes->globalDegreesOfFreedom;
52 }
53 myFirstDOF=Paso_Distribution_getFirstComponent(dof_distribution);
54 myLastDOF=Paso_Distribution_getLastComponent(dof_distribution);
55
56
57 mpiSize=mpi_info->size;
58 myRank=mpi_info->rank;
59
60 min_DOF=Finley_Util_getFlaggedMinInt(1,numNodes,globalDOFIndex,-1);
61 max_DOF=Finley_Util_getFlaggedMaxInt(1,numNodes,globalDOFIndex,-1);
62
63 if (max_DOF < min_DOF) {
64 min_DOF=myFirstDOF;
65 max_DOF=myLastDOF-1;
66 }
67
68 p_min=mpiSize;
69 p_max=-1;
70 if (max_DOF >= min_DOF) {
71 for (p=0; p<mpiSize; ++p) {
72 if (dof_distribution->first_component[p]<=min_DOF) p_min=p;
73 if (dof_distribution->first_component[p]<=max_DOF) p_max=p;
74 }
75 }
76
77 len_loc_dof=max_DOF-min_DOF+1;
78 if (! ((min_DOF<=myFirstDOF) && (myLastDOF-1<=max_DOF)) ) {
79 Finley_setError(SYSTEM_ERROR,"Local elements do not span local degrees of freedom.");
80 return;
81 }
82 rcv_len=TMPMEMALLOC(mpiSize,dim_t);
83 snd_len=TMPMEMALLOC(mpiSize,dim_t);
84 #ifdef PASO_MPI
85 mpi_requests=MEMALLOC(mpiSize*2,MPI_Request);
86 mpi_stati=MEMALLOC(mpiSize*2,MPI_Status);
87 #else
88 mpi_requests=MEMALLOC(mpiSize*2,int);
89 mpi_stati=MEMALLOC(mpiSize*2,int);
90 #endif
91 wanted_DOFs=TMPMEMALLOC(numNodes,index_t);
92 nodeMask=TMPMEMALLOC(numNodes,index_t);
93 neighbor=TMPMEMALLOC(mpiSize,Paso_MPI_rank);
94 shared=TMPMEMALLOC(numNodes*(p_max-p_min+1),index_t);
95 offsetInShared=TMPMEMALLOC(mpiSize+1,index_t);
96 locDOFMask=TMPMEMALLOC(len_loc_dof, index_t);
97 if (! ( Finley_checkPtr(neighbor) || Finley_checkPtr(shared) || Finley_checkPtr(offsetInShared) || Finley_checkPtr(locDOFMask) ||
98 Finley_checkPtr(nodeMask) || Finley_checkPtr(rcv_len) || Finley_checkPtr(snd_len) || Finley_checkPtr(mpi_requests) || Finley_checkPtr(mpi_stati) ||
99 Finley_checkPtr(mpi_stati) )) {
100
101 memset(rcv_len,0,sizeof(dim_t)*mpiSize);
102 #pragma omp parallel
103 {
104 #pragma omp for private(i) schedule(static)
105 for (i=0;i<len_loc_dof;++i) locDOFMask[i]=UNUSED;
106 #pragma omp for private(i) schedule(static)
107 for (i=0;i<numNodes;++i) nodeMask[i]=UNUSED;
108 #pragma omp for private(i,k) schedule(static)
109 for (i=0;i<numNodes;++i) {
110 k=globalDOFIndex[i];
111 if (k>-1) {
112 locDOFMask[k-min_DOF]=UNUSED-1;
113 #ifdef BOUNDS_CHECK
114 if ((k-min_DOF) >= len_loc_dof) { printf("BOUNDS_CHECK %s %d i=%d k=%d min_DOF=%d\n", __FILE__, __LINE__, i, k, min_DOF); exit(1); }
115 #endif
116 }
117 }
118
119 #pragma omp for private(i) schedule(static)
120 for (i=myFirstDOF-min_DOF;i<myLastDOF-min_DOF;++i) {
121 locDOFMask[i]=i-myFirstDOF+min_DOF;
122 #ifdef BOUNDS_CHECK
123 if (i < 0 || i >= len_loc_dof) { printf("BOUNDS_CHECK %s %d i=%d\n", __FILE__, __LINE__, i); exit(1); }
124 #endif
125 }
126 }
127
128 numNeighbors=0;
129 n=0;
130 lastn=n;
131 for (p=p_min;p<=p_max;++p) {
132 firstDOF=MAX(min_DOF,dof_distribution->first_component[p]);
133 lastDOF=MIN(max_DOF+1,dof_distribution->first_component[p+1]);
134 if (p != myRank) {
135 for (i=firstDOF-min_DOF;i<lastDOF-min_DOF;++i) {
136 #ifdef BOUNDS_CHECK
137 if (i < 0 || i >= len_loc_dof) { printf("BOUNDS_CHECK %s %d p=%d i=%d\n", __FILE__, __LINE__, p, i); exit(1); }
138 #endif
139 if (locDOFMask[i] == UNUSED-1) {
140 locDOFMask[i]=myLastDOF-myFirstDOF+n;
141 wanted_DOFs[n]=i+min_DOF;
142 ++n;
143 }
144 }
145 if (n>lastn) {
146 rcv_len[p]=n-lastn;
147 neighbor[numNeighbors]=p;
148 #ifdef BOUNDS_CHECK
149 if (numNeighbors < 0 || numNeighbors >= mpiSize+1) { printf("BOUNDS_CHECK %s %d p=%d numNeighbors=%d n=%d\n", __FILE__, __LINE__, p, numNeighbors, n); exit(1); }
150 #endif
151 offsetInShared[numNeighbors]=lastn;
152 numNeighbors++;
153 lastn=n;
154 }
155 }
156 }
157 #ifdef BOUNDS_CHECK
158 if (numNeighbors < 0 || numNeighbors >= mpiSize+1) { printf("BOUNDS_CHECK %s %d numNeighbors=%d\n", __FILE__, __LINE__, numNeighbors); exit(1); }
159 #endif
160 offsetInShared[numNeighbors]=lastn;
161
162 /* assign new DOF labels to nodes */
163 #pragma omp parallel for private(i) schedule(static)
164 for (i=0;i<numNodes;++i) {
165 k=globalDOFIndex[i];
166 if (k>-1) nodeMask[i]=locDOFMask[k-min_DOF];
167 }
168
169 /* now we can set the mapping from nodes to local DOFs */
170 this_mapping=Finley_NodeMapping_alloc(numNodes,nodeMask,UNUSED);
171 /* define how to get DOF values for controlled bu other processors */
172 #ifdef BOUNDS_CHECK
173 for (i=0;i<offsetInShared[numNeighbors];++i) {
174 if (i < 0 || i >= numNodes*(p_max-p_min+1)) { printf("BOUNDS_CHECK %s %d i=%d\n", __FILE__, __LINE__, i); exit(1); }
175 }
176 #endif
177 #pragma omp parallel for private(i) schedule(static)
178 for (i=0;i<offsetInShared[numNeighbors];++i) shared[i]=myLastDOF-myFirstDOF+i;
179
180 rcv_shcomp=Paso_SharedComponents_alloc(myLastDOF-myFirstDOF,numNeighbors,neighbor,shared,offsetInShared,1,0,mpi_info);
181
182 /*
183 * now we build the sender
184 */
185 #ifdef PASO_MPI
186 MPI_Alltoall(rcv_len,1,MPI_INT,snd_len,1,MPI_INT,mpi_info->comm);
187 #else
188 for (p=0;p<mpiSize;++p) snd_len[p]=rcv_len[p];
189 #endif
190 count=0;
191 for (p=0;p<rcv_shcomp->numNeighbors;p++) {
192 #ifdef PASO_MPI
193 MPI_Isend(&(wanted_DOFs[rcv_shcomp->offsetInShared[p]]), rcv_shcomp->offsetInShared[p+1]-rcv_shcomp->offsetInShared[p],
194 MPI_INT,rcv_shcomp->neighbor[p],mpi_info->msg_tag_counter+myRank,mpi_info->comm,&mpi_requests[count]);
195 #endif
196 count++;
197 }
198 n=0;
199 numNeighbors=0;
200 for (p=0;p<mpiSize;p++) {
201 if (snd_len[p] > 0) {
202 #ifdef PASO_MPI
203 MPI_Irecv(&(shared[n]),snd_len[p],
204 MPI_INT, p, mpi_info->msg_tag_counter+p, mpi_info->comm, &mpi_requests[count]);
205 #endif
206 count++;
207 neighbor[numNeighbors]=p;
208 offsetInShared[numNeighbors]=n;
209 numNeighbors++;
210 n+=snd_len[p];
211 }
212 }
213 mpi_info->msg_tag_counter+=mpi_info->size;
214 offsetInShared[numNeighbors]=n;
215 #ifdef PASO_MPI
216 MPI_Waitall(count,mpi_requests,mpi_stati);
217 #endif
218 /* map global ids to local id's */
219 #pragma omp parallel for private(i) schedule(static)
220 for (i=0;i<offsetInShared[numNeighbors];++i) {
221 shared[i]=locDOFMask[shared[i]-min_DOF];
222 }
223
224 snd_shcomp=Paso_SharedComponents_alloc(myLastDOF-myFirstDOF,numNeighbors,neighbor,shared,offsetInShared,1,0,dof_distribution->mpi_info);
225
226 if (Finley_noError()) this_connector=Paso_Connector_alloc(snd_shcomp,rcv_shcomp);
227 /* assign new DOF labels to nodes */
228 Paso_SharedComponents_free(rcv_shcomp);
229 Paso_SharedComponents_free(snd_shcomp);
230 }
231 TMPMEMFREE(rcv_len);
232 TMPMEMFREE(snd_len);
233 TMPMEMFREE(mpi_requests);
234 TMPMEMFREE(mpi_stati);
235 TMPMEMFREE(wanted_DOFs);
236 TMPMEMFREE(nodeMask);
237 TMPMEMFREE(neighbor);
238 TMPMEMFREE(shared);
239 TMPMEMFREE(offsetInShared);
240 TMPMEMFREE(locDOFMask);
241 if (Finley_noError()) {
242 if (use_reduced_elements) {
243 in->Nodes->reducedDegreesOfFreedomMapping=this_mapping;
244 in->Nodes->reducedDegreesOfFreedomConnector=this_connector;
245 } else {
246 in->Nodes->degreesOfFreedomMapping=this_mapping;
247 in->Nodes->degreesOfFreedomConnector=this_connector;
248 }
249 } else {
250 Finley_NodeMapping_free(this_mapping);
251 Paso_Connector_free(this_connector);
252
253 }
254 }
255
256 void Finley_Mesh_createMappings(Finley_Mesh* mesh, index_t* dof_distribution, index_t* node_distribution) {
257 int i;
258 index_t *maskReducedNodes=NULL, *indexReducedNodes=NULL;
259 dim_t numReducedNodes;
260
261 maskReducedNodes=TMPMEMALLOC(mesh->Nodes->numNodes,index_t);
262 indexReducedNodes=TMPMEMALLOC(mesh->Nodes->numNodes,index_t);
263
264 if (! ( Finley_checkPtr(maskReducedNodes) || Finley_checkPtr(indexReducedNodes) ) ) {
265 #pragma omp parallel for private(i) schedule(static)
266 for (i=0;i<mesh->Nodes->numNodes;++i) maskReducedNodes[i]=-1;
267 Finley_Mesh_markNodes(maskReducedNodes,0,mesh,TRUE);
268 numReducedNodes=Finley_Util_packMask(mesh->Nodes->numNodes,maskReducedNodes,indexReducedNodes);
269 if (Finley_noError()) Finley_Mesh_createNodeFileMappings(mesh,numReducedNodes,indexReducedNodes,dof_distribution, node_distribution);
270 }
271
272 TMPMEMFREE(maskReducedNodes);
273 TMPMEMFREE(indexReducedNodes);
274 }
275
276 void Finley_Mesh_createNodeFileMappings(Finley_Mesh* in, dim_t numReducedNodes, index_t* indexReducedNodes, index_t* dof_first_component, index_t* nodes_first_component) {
277
278
279 index_t myFirstDOF, myLastDOF, myFirstNode, myLastNode, *reduced_dof_first_component=NULL, *nodeMask=NULL,
280 *reduced_nodes_first_component=NULL, k,*maskMyReducedDOF=NULL, *indexMyReducedDOF=NULL, *maskMyReducedNodes=NULL, *indexMyReducedNodes=NULL;
281 dim_t myNumDOF, myNumNodes, myNumReducedNodes, myNumReducedDOF, globalNumReducedNodes, globalNumReducedDOF,i,mpiSize, minGlobalNodeIndex,maxGlobalNodeIndex;
282 Paso_MPI_rank myRank;
283
284 mpiSize=in->Nodes->MPIInfo->size;
285 myRank=in->Nodes->MPIInfo->rank;
286 /* mark the nodes used by the reduced mesh */
287
288 reduced_dof_first_component=TMPMEMALLOC(mpiSize+1,index_t);
289 reduced_nodes_first_component=TMPMEMALLOC(mpiSize+1,index_t);
290
291 if (! ( Finley_checkPtr(reduced_dof_first_component) || Finley_checkPtr(reduced_nodes_first_component) ) ) {
292
293 myFirstDOF=dof_first_component[myRank];
294 myLastDOF=dof_first_component[myRank+1];
295 myNumDOF=myLastDOF-myFirstDOF;
296
297 myFirstNode=nodes_first_component[myRank];
298 myLastNode=nodes_first_component[myRank+1];
299 myNumNodes=myLastNode-myFirstNode;
300
301 maskMyReducedDOF=TMPMEMALLOC(myNumDOF,index_t);
302 indexMyReducedDOF=TMPMEMALLOC(myNumDOF,index_t);
303 maskMyReducedNodes=TMPMEMALLOC(myNumNodes,index_t);
304 indexMyReducedNodes=TMPMEMALLOC(myNumNodes,index_t);
305
306 if (! ( Finley_checkPtr(maskMyReducedDOF) || Finley_checkPtr(indexMyReducedDOF) || Finley_checkPtr(maskMyReducedNodes) || Finley_checkPtr(indexMyReducedNodes) ) ) {
307
308 #pragma omp parallel private(i)
309 {
310 #pragma omp for schedule(static)
311 for (i=0;i<myNumNodes;++i) maskMyReducedNodes[i]=-1;
312 #pragma omp for schedule(static)
313 for (i=0;i<myNumDOF;++i) maskMyReducedDOF[i]=-1;
314 #pragma omp for private(k) schedule(static)
315 for (i=0;i<numReducedNodes;++i) {
316 k=in->Nodes->globalNodesIndex[indexReducedNodes[i]];
317 if ( (k>=myFirstNode) && (myLastNode>k) ) maskMyReducedNodes[k-myFirstNode]=i;
318 k=in->Nodes->globalDegreesOfFreedom[indexReducedNodes[i]];
319 if ( (k>=myFirstDOF) && (myLastDOF>k) ) {
320 maskMyReducedDOF[k-myFirstDOF]=i;
321 }
322 }
323 }
324 myNumReducedNodes=Finley_Util_packMask(myNumNodes,maskMyReducedNodes,indexMyReducedNodes);
325 myNumReducedDOF=Finley_Util_packMask(myNumDOF,maskMyReducedDOF,indexMyReducedDOF);
326
327 #ifdef PASO_MPI
328 MPI_Allgather(&myNumReducedNodes,1,MPI_INT,reduced_nodes_first_component,1,MPI_INT,in->Nodes->MPIInfo->comm);
329 MPI_Allgather(&myNumReducedDOF,1,MPI_INT,reduced_dof_first_component,1,MPI_INT,in->Nodes->MPIInfo->comm);
330 #else
331 reduced_nodes_first_component[0]=myNumReducedNodes;
332 reduced_dof_first_component[0]=myNumReducedDOF;
333 #endif
334 globalNumReducedNodes=0;
335 globalNumReducedDOF=0;
336 for (i=0;i<mpiSize;++i) {
337 k=reduced_nodes_first_component[i];
338 reduced_nodes_first_component[i]=globalNumReducedNodes;
339 globalNumReducedNodes+=k;
340
341 k=reduced_dof_first_component[i];
342 reduced_dof_first_component[i]=globalNumReducedDOF;
343 globalNumReducedDOF+=k;
344 }
345 reduced_nodes_first_component[mpiSize]=globalNumReducedNodes;
346 reduced_dof_first_component[mpiSize]=globalNumReducedDOF;
347 /* ==== distribution of Nodes ===============================*/
348 in->Nodes->nodesDistribution=Paso_Distribution_alloc(in->Nodes->MPIInfo,nodes_first_component,1,0);
349
350 /* ==== distribution of DOFs ===============================*/
351 in->Nodes->degreesOfFreedomDistribution=Paso_Distribution_alloc(in->Nodes->MPIInfo,dof_first_component,1,0);
352
353 /* ==== distribution of reduced Nodes ===============================*/
354 in->Nodes->reducedNodesDistribution=Paso_Distribution_alloc(in->Nodes->MPIInfo,reduced_nodes_first_component,1,0);
355
356 /* ==== distribution of reduced DOF ===============================*/
357 in->Nodes->reducedDegreesOfFreedomDistribution=Paso_Distribution_alloc(in->Nodes->MPIInfo,reduced_dof_first_component,1,0);
358 }
359 TMPMEMFREE(maskMyReducedDOF);
360 TMPMEMFREE(indexMyReducedDOF);
361 TMPMEMFREE(maskMyReducedNodes);
362 TMPMEMFREE(indexMyReducedNodes);
363 }
364 TMPMEMFREE(reduced_dof_first_component);
365 TMPMEMFREE(reduced_nodes_first_component);
366
367 nodeMask=TMPMEMALLOC(in->Nodes->numNodes,index_t);
368 if (! Finley_checkPtr(nodeMask) && Finley_noError()) {
369
370 /* ==== nodes mapping which is a dummy structure ======== */
371 #pragma omp parallel for private(i) schedule(static)
372 for (i=0;i<in->Nodes->numNodes;++i) nodeMask[i]=i;
373 in->Nodes->nodesMapping=Finley_NodeMapping_alloc(in->Nodes->numNodes,nodeMask,UNUSED);
374
375 /* ==== mapping between nodes and reduced nodes ========== */
376 #pragma omp parallel for private(i) schedule(static)
377 for (i=0;i<in->Nodes->numNodes;++i) nodeMask[i]=UNUSED;
378 #pragma omp parallel for private(i) schedule(static)
379 for (i=0;i<numReducedNodes;++i) nodeMask[indexReducedNodes[i]]=i;
380 in->Nodes->reducedNodesMapping=Finley_NodeMapping_alloc(in->Nodes->numNodes,nodeMask,UNUSED);
381
382 }
383 TMPMEMFREE(nodeMask);
384 /* ==== mapping between nodes and DOFs + DOF connector ========== */
385 if ( Finley_noError()) Mesh_createDOFMappingAndCoupling(in,FALSE);
386 /* ==== mapping between nodes and reduced DOFs + reduced DOF connector ========== */
387 if ( Finley_noError()) Mesh_createDOFMappingAndCoupling(in,TRUE);
388
389 /* get the Ids for DOFs and reduced nodes */
390 if (Finley_noError()) {
391 #pragma omp parallel private(i)
392 {
393 #pragma omp for
394 for (i=0;i<in->Nodes->reducedNodesMapping->numTargets;++i) in->Nodes->reducedNodesId[i]=in->Nodes->Id[in->Nodes->reducedNodesMapping->map[i]];
395 #pragma omp for
396 for (i=0;i<in->Nodes->degreesOfFreedomMapping->numTargets;++i) in->Nodes->degreesOfFreedomId[i]=in->Nodes->Id[in->Nodes->degreesOfFreedomMapping->map[i]];
397 #pragma omp for
398 for (i=0;i<in->Nodes->reducedDegreesOfFreedomMapping->numTargets;++i) in->Nodes->reducedDegreesOfFreedomId[i]=in->Nodes->Id[in->Nodes->reducedDegreesOfFreedomMapping->map[i]];
399 }
400 } else {
401 Finley_NodeMapping_free(in->Nodes->nodesMapping);
402 Finley_NodeMapping_free(in->Nodes->reducedNodesMapping);
403 Finley_NodeMapping_free(in->Nodes->degreesOfFreedomMapping);
404 Finley_NodeMapping_free(in->Nodes->reducedDegreesOfFreedomMapping);
405 Paso_Distribution_free(in->Nodes->nodesDistribution);
406 Paso_Distribution_free(in->Nodes->reducedNodesDistribution);
407 Paso_Distribution_free(in->Nodes->degreesOfFreedomDistribution);
408 Paso_Distribution_free(in->Nodes->reducedDegreesOfFreedomDistribution);
409 Paso_Connector_free(in->Nodes->degreesOfFreedomConnector);
410 Paso_Connector_free(in->Nodes->reducedDegreesOfFreedomConnector);
411 in->Nodes->nodesMapping=NULL;
412 in->Nodes->reducedNodesMapping=NULL;
413 in->Nodes->degreesOfFreedomMapping=NULL;
414 in->Nodes->reducedDegreesOfFreedomMapping=NULL;
415 in->Nodes->nodesDistribution=NULL;
416 in->Nodes->reducedNodesDistribution=NULL;
417 in->Nodes->degreesOfFreedomDistribution=NULL;
418 in->Nodes->reducedDegreesOfFreedomDistribution=NULL;
419 in->Nodes->degreesOfFreedomConnector=NULL;
420 in->Nodes->reducedDegreesOfFreedomConnector=NULL;
421 }
422 }
423

  ViewVC Help
Powered by ViewVC 1.1.26