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