/[escript]/trunk/finley/src/CPPAdapter/MeshAdapter.cpp
ViewVC logotype

Contents of /trunk/finley/src/CPPAdapter/MeshAdapter.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1347 - (show annotations)
Fri Nov 16 05:37:07 2007 UTC (11 years, 7 months ago) by ksteube
File size: 82429 byte(s)
Completed mesh.dump(file) and mesh=LoadMesh(file) by adding TagMap and
implementing MPI parallelism.
Now allocating ElementFile for ContactElements even if there are none.
Removed file Mesh_dump.c since dump/loadMesh are in CPPAdapter/MeshAdapter*.cpp.

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 #include "MeshAdapter.h"
17 #include "escript/Data.h"
18 #include "escript/DataFactory.h"
19 #ifdef USE_NETCDF
20 #include <netcdfcpp.h>
21 #endif
22 extern "C" {
23 #include "escript/blocktimer.h"
24 }
25 #include <vector>
26
27 using namespace std;
28 using namespace escript;
29
30 namespace finley {
31
32 //
33 // define the static constants
34 MeshAdapter::FunctionSpaceNamesMapType MeshAdapter::m_functionSpaceTypeNames;
35 const int MeshAdapter::DegreesOfFreedom=FINLEY_DEGREES_OF_FREEDOM;
36 const int MeshAdapter::ReducedDegreesOfFreedom=FINLEY_REDUCED_DEGREES_OF_FREEDOM;
37 const int MeshAdapter::Nodes=FINLEY_NODES;
38 const int MeshAdapter::ReducedNodes=FINLEY_REDUCED_NODES;
39 const int MeshAdapter::Elements=FINLEY_ELEMENTS;
40 const int MeshAdapter::ReducedElements=FINLEY_REDUCED_ELEMENTS;
41 const int MeshAdapter::FaceElements=FINLEY_FACE_ELEMENTS;
42 const int MeshAdapter::ReducedFaceElements=FINLEY_REDUCED_FACE_ELEMENTS;
43 const int MeshAdapter::Points=FINLEY_POINTS;
44 const int MeshAdapter::ContactElementsZero=FINLEY_CONTACT_ELEMENTS_1;
45 const int MeshAdapter::ReducedContactElementsZero=FINLEY_REDUCED_CONTACT_ELEMENTS_1;
46 const int MeshAdapter::ContactElementsOne=FINLEY_CONTACT_ELEMENTS_2;
47 const int MeshAdapter::ReducedContactElementsOne=FINLEY_REDUCED_CONTACT_ELEMENTS_2;
48
49 MeshAdapter::MeshAdapter(Finley_Mesh* finleyMesh)
50 {
51 setFunctionSpaceTypeNames();
52 //
53 // need to use a null_deleter as Finley_Mesh_free deletes the pointer
54 // for us.
55 m_finleyMesh.reset(finleyMesh,null_deleter());
56 }
57
58 //
59 // The copy constructor should just increment the use count
60 MeshAdapter::MeshAdapter(const MeshAdapter& in):
61 m_finleyMesh(in.m_finleyMesh)
62 {
63 setFunctionSpaceTypeNames();
64 }
65
66 MeshAdapter::~MeshAdapter()
67 {
68 //
69 // I hope the case for the pointer being zero has been taken care of.
70 // cout << "In MeshAdapter destructor." << endl;
71 if (m_finleyMesh.unique()) {
72 Finley_Mesh_free(m_finleyMesh.get());
73 }
74 }
75
76 int MeshAdapter::getMPISize() const
77 {
78 return m_finleyMesh.get()->MPIInfo->size;
79 }
80 int MeshAdapter::getMPIRank() const
81 {
82 return m_finleyMesh.get()->MPIInfo->rank;
83 }
84
85
86 Finley_Mesh* MeshAdapter::getFinley_Mesh() const {
87 return m_finleyMesh.get();
88 }
89
90 void MeshAdapter::write(const std::string& fileName) const
91 {
92 char *fName = (fileName.size()+1>0) ? TMPMEMALLOC(fileName.size()+1,char) : (char*)NULL;
93 strcpy(fName,fileName.c_str());
94 Finley_Mesh_write(m_finleyMesh.get(),fName);
95 checkFinleyError();
96 TMPMEMFREE(fName);
97 }
98
99 void MeshAdapter::Print_Mesh_Info(const bool full=false) const
100 {
101 Finley_PrintMesh_Info(m_finleyMesh.get(), full);
102 }
103
104 void MeshAdapter::dump(const std::string& fileName) const
105 {
106 #ifdef USE_NETCDF
107 const NcDim* ncdims[12];
108 NcVar *ids, *data;
109 int *int_ptr;
110 Finley_Mesh *mesh = m_finleyMesh.get();
111 Finley_TagMap* tag_map;
112 int num_Tags = 0;
113 int mpi_size = mesh->MPIInfo->size;
114 int mpi_rank = mesh->MPIInfo->rank;
115 int numDim = mesh->Nodes->numDim;
116 int numNodes = mesh->Nodes->numNodes;
117 int num_Elements = mesh->Elements->numElements;
118 int num_FaceElements = mesh->FaceElements->numElements;
119 int num_ContactElements = mesh->ContactElements->numElements;
120 int num_Points = mesh->Points->numElements;
121 int num_Elements_numNodes = mesh->Elements->numNodes;
122 int num_FaceElements_numNodes = mesh->FaceElements->numNodes;
123 int num_ContactElements_numNodes = mesh->ContactElements->numNodes;
124 char *newFileName = Paso_MPI_appendRankToFileName(strdup(fileName.c_str()), mpi_size, mpi_rank);
125
126 /* Figure out how much storage is required for tags */
127 tag_map = mesh->TagMap;
128 if (tag_map) {
129 while (tag_map) {
130 num_Tags++;
131 tag_map=tag_map->next;
132 }
133 }
134
135 // NetCDF error handler
136 NcError err(NcError::verbose_nonfatal);
137 // Create the file.
138 NcFile dataFile(newFileName, NcFile::Replace);
139 // check if writing was successful
140 if (!dataFile.is_valid())
141 throw DataException("Error - MeshAdapter::dump: opening of NetCDF file for output failed: " + *newFileName);
142
143 // Define dimensions (num_Elements and dim_Elements are identical, dim_Elements only appears if > 0)
144 if (! (ncdims[0] = dataFile.add_dim("numNodes", numNodes)) )
145 throw DataException("Error - MeshAdapter::dump: appending dimension numNodes to netCDF file failed: " + *newFileName);
146 if (! (ncdims[1] = dataFile.add_dim("numDim", numDim)) )
147 throw DataException("Error - MeshAdapter::dump: appending dimension numDim to netCDF file failed: " + *newFileName);
148 if (! (ncdims[2] = dataFile.add_dim("mpi_size_plus_1", mpi_size+1)) )
149 throw DataException("Error - MeshAdapter::dump: appending dimension mpi_size to netCDF file failed: " + *newFileName);
150 if (num_Elements>0)
151 if (! (ncdims[3] = dataFile.add_dim("dim_Elements", num_Elements)) )
152 throw DataException("Error - MeshAdapter::dump: appending dimension dim_Elements to netCDF file failed: " + *newFileName);
153 if (num_FaceElements>0)
154 if (! (ncdims[4] = dataFile.add_dim("dim_FaceElements", num_FaceElements)) )
155 throw DataException("Error - MeshAdapter::dump: appending dimension dim_FaceElements to netCDF file failed: " + *newFileName);
156 if (num_ContactElements>0)
157 if (! (ncdims[5] = dataFile.add_dim("dim_ContactElements", num_ContactElements)) )
158 throw DataException("Error - MeshAdapter::dump: appending dimension dim_ContactElements to netCDF file failed: " + *newFileName);
159 if (num_Points>0)
160 if (! (ncdims[6] = dataFile.add_dim("dim_Points", num_Points)) )
161 throw DataException("Error - MeshAdapter::dump: appending dimension dim_Points to netCDF file failed: " + *newFileName);
162 if (num_Elements>0)
163 if (! (ncdims[7] = dataFile.add_dim("dim_Elements_Nodes", num_Elements_numNodes)) )
164 throw DataException("Error - MeshAdapter::dump: appending dimension dim_Elements_Nodes to netCDF file failed: " + *newFileName);
165 if (num_FaceElements>0)
166 if (! (ncdims[8] = dataFile.add_dim("dim_FaceElements_numNodes", num_FaceElements_numNodes)) )
167 throw DataException("Error - MeshAdapter::dump: appending dimension dim_FaceElements_numNodes to netCDF file failed: " + *newFileName);
168 if (num_ContactElements>0)
169 if (! (ncdims[9] = dataFile.add_dim("dim_ContactElements_numNodes", num_ContactElements_numNodes)) )
170 throw DataException("Error - MeshAdapter::dump: appending dimension dim_ContactElements_numNodes to netCDF file failed: " + *newFileName);
171 if (num_Tags>0)
172 if (! (ncdims[10] = dataFile.add_dim("dim_Tags", num_Tags)) )
173 throw DataException("Error - MeshAdapter::dump: appending dimension dim_Tags to netCDF file failed: " + *newFileName);
174
175 // Attributes: MPI size, MPI rank, Name, order, reduced_order
176 if (!dataFile.add_att("mpi_size", mpi_size) )
177 throw DataException("Error - MeshAdapter::dump: appending mpi_size to NetCDF file failed: " + *newFileName);
178 if (!dataFile.add_att("mpi_rank", mpi_rank) )
179 throw DataException("Error - MeshAdapter::dump: appending mpi_rank to NetCDF file failed: " + *newFileName);
180 if (!dataFile.add_att("Name",mesh->Name) )
181 throw DataException("Error - MeshAdapter::dump: appending Name to NetCDF file failed: " + *newFileName);
182 if (!dataFile.add_att("numDim",numDim) )
183 throw DataException("Error - MeshAdapter::dump: appending order to NetCDF file failed: " + *newFileName);
184 if (!dataFile.add_att("order",mesh->order) )
185 throw DataException("Error - MeshAdapter::dump: appending order to NetCDF file failed: " + *newFileName);
186 if (!dataFile.add_att("reduced_order",mesh->reduced_order) )
187 throw DataException("Error - MeshAdapter::dump: appending reduced_order to NetCDF file failed: " + *newFileName);
188 if (!dataFile.add_att("numNodes",numNodes) )
189 throw DataException("Error - MeshAdapter::dump: appending numNodes to NetCDF file failed: " + *newFileName);
190 if (!dataFile.add_att("num_Elements",num_Elements) )
191 throw DataException("Error - MeshAdapter::dump: appending num_Elements to NetCDF file failed: " + *newFileName);
192 if (!dataFile.add_att("num_FaceElements",num_FaceElements) )
193 throw DataException("Error - MeshAdapter::dump: appending num_FaceElements to NetCDF file failed: " + *newFileName);
194 if (!dataFile.add_att("num_ContactElements",num_ContactElements) )
195 throw DataException("Error - MeshAdapter::dump: appending num_ContactElements to NetCDF file failed: " + *newFileName);
196 if (!dataFile.add_att("num_Points",num_Points) )
197 throw DataException("Error - MeshAdapter::dump: appending num_Points to NetCDF file failed: " + *newFileName);
198 if (!dataFile.add_att("num_Elements_numNodes",num_Elements_numNodes) )
199 throw DataException("Error - MeshAdapter::dump: appending num_Elements_numNodes to NetCDF file failed: " + *newFileName);
200 if (!dataFile.add_att("num_FaceElements_numNodes",num_FaceElements_numNodes) )
201 throw DataException("Error - MeshAdapter::dump: appending num_FaceElements_numNodes to NetCDF file failed: " + *newFileName);
202 if (!dataFile.add_att("num_ContactElements_numNodes",num_ContactElements_numNodes) )
203 throw DataException("Error - MeshAdapter::dump: appending num_ContactElements_numNodes to NetCDF file failed: " + *newFileName);
204 if (!dataFile.add_att("Elements_TypeId", mesh->Elements->ReferenceElement->Type->TypeId) )
205 throw DataException("Error - MeshAdapter::dump: appending Elements_TypeId to NetCDF file failed: " + *newFileName);
206 if (!dataFile.add_att("FaceElements_TypeId", mesh->FaceElements->ReferenceElement->Type->TypeId) )
207 throw DataException("Error - MeshAdapter::dump: appending FaceElements_TypeId to NetCDF file failed: " + *newFileName);
208 if (!dataFile.add_att("ContactElements_TypeId", mesh->ContactElements->ReferenceElement->Type->TypeId) )
209 throw DataException("Error - MeshAdapter::dump: appending ContactElements_TypeId to NetCDF file failed: " + *newFileName);
210 if (!dataFile.add_att("Points_TypeId", mesh->Points->ReferenceElement->Type->TypeId) )
211 throw DataException("Error - MeshAdapter::dump: appending Points_TypeId to NetCDF file failed: " + *newFileName);
212 if (!dataFile.add_att("num_Tags", num_Tags) )
213 throw DataException("Error - MeshAdapter::dump: appending num_Tags to NetCDF file failed: " + *newFileName);
214
215 // // // // // Nodes // // // // //
216
217 // Only write nodes if non-empty because NetCDF doesn't like empty arrays (it treats them as NC_UNLIMITED)
218 if (numNodes>0) {
219
220 // Nodes Id
221 if (! ( ids = dataFile.add_var("Nodes_Id", ncInt, ncdims[0])) )
222 throw DataException("Error - MeshAdapter::dump: appending Nodes_Id to netCDF file failed: " + *newFileName);
223 int_ptr = &mesh->Nodes->Id[0];
224 if (! (ids->put(int_ptr, numNodes)) )
225 throw DataException("Error - MeshAdapter::dump: copy Nodes_Id to netCDF buffer failed: " + *newFileName);
226
227 // Nodes Tag
228 if (! ( ids = dataFile.add_var("Nodes_Tag", ncInt, ncdims[0])) )
229 throw DataException("Error - MeshAdapter::dump: appending Nodes_Tag to netCDF file failed: " + *newFileName);
230 int_ptr = &mesh->Nodes->Tag[0];
231 if (! (ids->put(int_ptr, numNodes)) )
232 throw DataException("Error - MeshAdapter::dump: copy Nodes_Tag to netCDF buffer failed: " + *newFileName);
233
234 // Nodes gDOF
235 if (! ( ids = dataFile.add_var("Nodes_gDOF", ncInt, ncdims[0])) )
236 throw DataException("Error - MeshAdapter::dump: appending Nodes_gDOF to netCDF file failed: " + *newFileName);
237 int_ptr = &mesh->Nodes->globalDegreesOfFreedom[0];
238 if (! (ids->put(int_ptr, numNodes)) )
239 throw DataException("Error - MeshAdapter::dump: copy Nodes_gDOF to netCDF buffer failed: " + *newFileName);
240
241 // Nodes global node index
242 if (! ( ids = dataFile.add_var("Nodes_gNI", ncInt, ncdims[0])) )
243 throw DataException("Error - MeshAdapter::dump: appending Nodes_gNI to netCDF file failed: " + *newFileName);
244 int_ptr = &mesh->Nodes->globalNodesIndex[0];
245 if (! (ids->put(int_ptr, numNodes)) )
246 throw DataException("Error - MeshAdapter::dump: copy Nodes_gNI to netCDF buffer failed: " + *newFileName);
247
248 // Nodes grDof
249 if (! ( ids = dataFile.add_var("Nodes_grDfI", ncInt, ncdims[0])) )
250 throw DataException("Error - MeshAdapter::dump: appending Nodes_grDfI to netCDF file failed: " + *newFileName);
251 int_ptr = &mesh->Nodes->globalReducedDOFIndex[0];
252 if (! (ids->put(int_ptr, numNodes)) )
253 throw DataException("Error - MeshAdapter::dump: copy Nodes_grDfI to netCDF buffer failed: " + *newFileName);
254
255 // Nodes grNI
256 if (! ( ids = dataFile.add_var("Nodes_grNI", ncInt, ncdims[0])) )
257 throw DataException("Error - MeshAdapter::dump: appending Nodes_grNI to netCDF file failed: " + *newFileName);
258 int_ptr = &mesh->Nodes->globalReducedNodesIndex[0];
259 if (! (ids->put(int_ptr, numNodes)) )
260 throw DataException("Error - MeshAdapter::dump: copy Nodes_grNI to netCDF buffer failed: " + *newFileName);
261
262 // Nodes Coordinates
263 if (! ( ids = dataFile.add_var("Nodes_Coordinates", ncDouble, ncdims[0], ncdims[1]) ) )
264 throw DataException("Error - MeshAdapter::dump: appending Nodes_Coordinates to netCDF file failed: " + *newFileName);
265 if (! (ids->put(&(mesh->Nodes->Coordinates[INDEX2(0,0,numDim)]), numNodes, numDim)) )
266 throw DataException("Error - MeshAdapter::dump: copy Nodes_Coordinates to netCDF buffer failed: " + *newFileName);
267
268 // Nodes degreesOfFreedomDistribution
269 if (! ( ids = dataFile.add_var("Nodes_DofDistribution", ncInt, ncdims[2])) )
270 throw DataException("Error - MeshAdapter::dump: appending Nodes_DofDistribution to netCDF file failed: " + *newFileName);
271 int_ptr = &mesh->Nodes->degreesOfFreedomDistribution->first_component[0];
272 if (! (ids->put(int_ptr, mpi_size+1)) )
273 throw DataException("Error - MeshAdapter::dump: copy Nodes_DofDistribution to netCDF buffer failed: " + *newFileName);
274
275 }
276
277 // // // // // Elements // // // // //
278
279 if (num_Elements>0) {
280
281 // Temp storage to gather node IDs
282 int *Elements_Nodes = TMPMEMALLOC(num_Elements*num_Elements_numNodes,int);
283
284 // Elements_Id
285 if (! ( ids = dataFile.add_var("Elements_Id", ncInt, ncdims[3])) )
286 throw DataException("Error - MeshAdapter::dump: appending Elements_Id to netCDF file failed: " + *newFileName);
287 int_ptr = &mesh->Elements->Id[0];
288 if (! (ids->put(int_ptr, num_Elements)) )
289 throw DataException("Error - MeshAdapter::dump: copy Elements_Id to netCDF buffer failed: " + *newFileName);
290
291 // Elements_Tag
292 if (! ( ids = dataFile.add_var("Elements_Tag", ncInt, ncdims[3])) )
293 throw DataException("Error - MeshAdapter::dump: appending Elements_Tag to netCDF file failed: " + *newFileName);
294 int_ptr = &mesh->Elements->Tag[0];
295 if (! (ids->put(int_ptr, num_Elements)) )
296 throw DataException("Error - MeshAdapter::dump: copy Elements_Tag to netCDF buffer failed: " + *newFileName);
297
298 // Elements_Owner
299 if (! ( ids = dataFile.add_var("Elements_Owner", ncInt, ncdims[3])) )
300 throw DataException("Error - MeshAdapter::dump: appending Elements_Owner to netCDF file failed: " + *newFileName);
301 int_ptr = &mesh->Elements->Owner[0];
302 if (! (ids->put(int_ptr, num_Elements)) )
303 throw DataException("Error - MeshAdapter::dump: copy Elements_Owner to netCDF buffer failed: " + *newFileName);
304
305 // Elements_Color
306 if (! ( ids = dataFile.add_var("Elements_Color", ncInt, ncdims[3])) )
307 throw DataException("Error - MeshAdapter::dump: appending Elements_Color to netCDF file failed: " + *newFileName);
308 int_ptr = &mesh->Elements->Color[0];
309 if (! (ids->put(int_ptr, num_Elements)) )
310 throw DataException("Error - MeshAdapter::dump: copy Elements_Color to netCDF buffer failed: " + *newFileName);
311
312 // Elements_Nodes
313 for (int i=0; i<num_Elements; i++)
314 for (int j=0; j<num_Elements_numNodes; j++)
315 Elements_Nodes[INDEX2(j,i,num_Elements_numNodes)] = mesh->Nodes->Id[mesh->Elements->Nodes[INDEX2(j,i,num_Elements_numNodes)]];
316 if (! ( ids = dataFile.add_var("Elements_Nodes", ncInt, ncdims[3], ncdims[7]) ) )
317 throw DataException("Error - MeshAdapter::dump: appending Elements_Nodes to netCDF file failed: " + *newFileName);
318 if (! (ids->put(&(Elements_Nodes[0]), num_Elements, num_Elements_numNodes)) )
319 throw DataException("Error - MeshAdapter::dump: copy Elements_Nodes to netCDF buffer failed: " + *newFileName);
320
321 TMPMEMFREE(Elements_Nodes);
322
323 }
324
325 // // // // // Face_Elements // // // // //
326
327 if (num_FaceElements>0) {
328
329 // Temp storage to gather node IDs
330 int *FaceElements_Nodes = TMPMEMALLOC(num_FaceElements*num_FaceElements_numNodes,int);
331
332 // FaceElements_Id
333 if (! ( ids = dataFile.add_var("FaceElements_Id", ncInt, ncdims[4])) )
334 throw DataException("Error - MeshAdapter::dump: appending FaceElements_Id to netCDF file failed: " + *newFileName);
335 int_ptr = &mesh->FaceElements->Id[0];
336 if (! (ids->put(int_ptr, num_FaceElements)) )
337 throw DataException("Error - MeshAdapter::dump: copy FaceElements_Id to netCDF buffer failed: " + *newFileName);
338
339 // FaceElements_Tag
340 if (! ( ids = dataFile.add_var("FaceElements_Tag", ncInt, ncdims[4])) )
341 throw DataException("Error - MeshAdapter::dump: appending FaceElements_Tag to netCDF file failed: " + *newFileName);
342 int_ptr = &mesh->FaceElements->Tag[0];
343 if (! (ids->put(int_ptr, num_FaceElements)) )
344 throw DataException("Error - MeshAdapter::dump: copy FaceElements_Tag to netCDF buffer failed: " + *newFileName);
345
346 // FaceElements_Owner
347 if (! ( ids = dataFile.add_var("FaceElements_Owner", ncInt, ncdims[4])) )
348 throw DataException("Error - MeshAdapter::dump: appending FaceElements_Owner to netCDF file failed: " + *newFileName);
349 int_ptr = &mesh->FaceElements->Owner[0];
350 if (! (ids->put(int_ptr, num_FaceElements)) )
351 throw DataException("Error - MeshAdapter::dump: copy FaceElements_Owner to netCDF buffer failed: " + *newFileName);
352
353 // FaceElements_Color
354 if (! ( ids = dataFile.add_var("FaceElements_Color", ncInt, ncdims[4])) )
355 throw DataException("Error - MeshAdapter::dump: appending FaceElements_Color to netCDF file failed: " + *newFileName);
356 int_ptr = &mesh->FaceElements->Color[0];
357 if (! (ids->put(int_ptr, num_FaceElements)) )
358 throw DataException("Error - MeshAdapter::dump: copy FaceElements_Color to netCDF buffer failed: " + *newFileName);
359
360 // FaceElements_Nodes
361 for (int i=0; i<num_FaceElements; i++)
362 for (int j=0; j<num_FaceElements_numNodes; j++)
363 FaceElements_Nodes[INDEX2(j,i,num_FaceElements_numNodes)] = mesh->Nodes->Id[mesh->FaceElements->Nodes[INDEX2(j,i,num_FaceElements_numNodes)]];
364 if (! ( ids = dataFile.add_var("FaceElements_Nodes", ncInt, ncdims[4], ncdims[8]) ) )
365 throw DataException("Error - MeshAdapter::dump: appending FaceElements_Nodes to netCDF file failed: " + *newFileName);
366 if (! (ids->put(&(FaceElements_Nodes[0]), num_FaceElements, num_FaceElements_numNodes)) )
367 throw DataException("Error - MeshAdapter::dump: copy FaceElements_Nodes to netCDF buffer failed: " + *newFileName);
368
369 TMPMEMFREE(FaceElements_Nodes);
370
371 }
372
373 // // // // // Contact_Elements // // // // //
374
375 if (num_ContactElements>0) {
376
377 // Temp storage to gather node IDs
378 int *ContactElements_Nodes = TMPMEMALLOC(num_ContactElements*num_ContactElements_numNodes,int);
379
380 // ContactElements_Id
381 if (! ( ids = dataFile.add_var("ContactElements_Id", ncInt, ncdims[5])) )
382 throw DataException("Error - MeshAdapter::dump: appending ContactElements_Id to netCDF file failed: " + *newFileName);
383 int_ptr = &mesh->ContactElements->Id[0];
384 if (! (ids->put(int_ptr, num_ContactElements)) )
385 throw DataException("Error - MeshAdapter::dump: copy ContactElements_Id to netCDF buffer failed: " + *newFileName);
386
387 // ContactElements_Tag
388 if (! ( ids = dataFile.add_var("ContactElements_Tag", ncInt, ncdims[5])) )
389 throw DataException("Error - MeshAdapter::dump: appending ContactElements_Tag to netCDF file failed: " + *newFileName);
390 int_ptr = &mesh->ContactElements->Tag[0];
391 if (! (ids->put(int_ptr, num_ContactElements)) )
392 throw DataException("Error - MeshAdapter::dump: copy ContactElements_Tag to netCDF buffer failed: " + *newFileName);
393
394 // ContactElements_Owner
395 if (! ( ids = dataFile.add_var("ContactElements_Owner", ncInt, ncdims[5])) )
396 throw DataException("Error - MeshAdapter::dump: appending ContactElements_Owner to netCDF file failed: " + *newFileName);
397 int_ptr = &mesh->ContactElements->Owner[0];
398 if (! (ids->put(int_ptr, num_ContactElements)) )
399 throw DataException("Error - MeshAdapter::dump: copy ContactElements_Owner to netCDF buffer failed: " + *newFileName);
400
401 // ContactElements_Color
402 if (! ( ids = dataFile.add_var("ContactElements_Color", ncInt, ncdims[5])) )
403 throw DataException("Error - MeshAdapter::dump: appending ContactElements_Color to netCDF file failed: " + *newFileName);
404 int_ptr = &mesh->ContactElements->Color[0];
405 if (! (ids->put(int_ptr, num_ContactElements)) )
406 throw DataException("Error - MeshAdapter::dump: copy ContactElements_Color to netCDF buffer failed: " + *newFileName);
407
408 // ContactElements_Nodes
409 for (int i=0; i<num_ContactElements; i++)
410 for (int j=0; j<num_ContactElements_numNodes; j++)
411 ContactElements_Nodes[INDEX2(j,i,num_ContactElements_numNodes)] = mesh->Nodes->Id[mesh->ContactElements->Nodes[INDEX2(j,i,num_ContactElements_numNodes)]];
412 if (! ( ids = dataFile.add_var("ContactElements_Nodes", ncInt, ncdims[5], ncdims[9]) ) )
413 throw DataException("Error - MeshAdapter::dump: appending ContactElements_Nodes to netCDF file failed: " + *newFileName);
414 if (! (ids->put(&(ContactElements_Nodes[0]), num_ContactElements, num_ContactElements_numNodes)) )
415 throw DataException("Error - MeshAdapter::dump: copy ContactElements_Nodes to netCDF buffer failed: " + *newFileName);
416
417 TMPMEMFREE(ContactElements_Nodes);
418
419 }
420
421 // // // // // Points // // // // //
422
423 if (num_Points>0) {
424
425 fprintf(stderr, "\n\n\nWARNING: MeshAdapter::dump has not been tested with Point elements\n\n\n");
426
427 // Temp storage to gather node IDs
428 int *Points_Nodes = TMPMEMALLOC(num_Points,int);
429
430 // Points_Id
431 if (! ( ids = dataFile.add_var("Points_Id", ncInt, ncdims[6])) )
432 throw DataException("Error - MeshAdapter::dump: appending Points_Id to netCDF file failed: " + *newFileName);
433 int_ptr = &mesh->Points->Id[0];
434 if (! (ids->put(int_ptr, num_Points)) )
435 throw DataException("Error - MeshAdapter::dump: copy Points_Id to netCDF buffer failed: " + *newFileName);
436
437 // Points_Tag
438 if (! ( ids = dataFile.add_var("Points_Tag", ncInt, ncdims[6])) )
439 throw DataException("Error - MeshAdapter::dump: appending Points_Tag to netCDF file failed: " + *newFileName);
440 int_ptr = &mesh->Points->Tag[0];
441 if (! (ids->put(int_ptr, num_Points)) )
442 throw DataException("Error - MeshAdapter::dump: copy Points_Tag to netCDF buffer failed: " + *newFileName);
443
444 // Points_Owner
445 if (! ( ids = dataFile.add_var("Points_Owner", ncInt, ncdims[6])) )
446 throw DataException("Error - MeshAdapter::dump: appending Points_Owner to netCDF file failed: " + *newFileName);
447 int_ptr = &mesh->Points->Owner[0];
448 if (! (ids->put(int_ptr, num_Points)) )
449 throw DataException("Error - MeshAdapter::dump: copy Points_Owner to netCDF buffer failed: " + *newFileName);
450
451 // Points_Color
452 if (! ( ids = dataFile.add_var("Points_Color", ncInt, ncdims[6])) )
453 throw DataException("Error - MeshAdapter::dump: appending Points_Color to netCDF file failed: " + *newFileName);
454 int_ptr = &mesh->Points->Color[0];
455 if (! (ids->put(int_ptr, num_Points)) )
456 throw DataException("Error - MeshAdapter::dump: copy Points_Color to netCDF buffer failed: " + *newFileName);
457
458 // Points_Nodes
459 // mesh->Nodes->Id[mesh->Points->Nodes[INDEX2(0,i,1)]]
460 for (int i=0; i<num_Points; i++)
461 Points_Nodes[i] = mesh->Nodes->Id[mesh->Points->Nodes[INDEX2(0,i,1)]];
462 if (! ( ids = dataFile.add_var("Points_Nodes", ncInt, ncdims[6]) ) )
463 throw DataException("Error - MeshAdapter::dump: appending Points_Nodes to netCDF file failed: " + *newFileName);
464 if (! (ids->put(&(Points_Nodes[0]), num_Points)) )
465 throw DataException("Error - MeshAdapter::dump: copy Points_Nodes to netCDF buffer failed: " + *newFileName);
466
467 TMPMEMFREE(Points_Nodes);
468
469 }
470
471 // // // // // TagMap // // // // //
472
473 if (num_Tags>0) {
474
475 // Temp storage to gather node IDs
476 int *Tags_keys = TMPMEMALLOC(num_Tags, int);
477 char name_temp[4096];
478
479 /* Copy tag data into temp arrays */
480 tag_map = mesh->TagMap;
481 if (tag_map) {
482 int i = 0;
483 while (tag_map) {
484 Tags_keys[i++] = tag_map->tag_key;
485 tag_map=tag_map->next;
486 }
487 }
488
489 // Tags_keys
490 if (! ( ids = dataFile.add_var("Tags_keys", ncInt, ncdims[10])) )
491 throw DataException("Error - MeshAdapter::dump: appending Tags_keys to netCDF file failed: " + *newFileName);
492 int_ptr = &Tags_keys[0];
493 if (! (ids->put(int_ptr, num_Tags)) )
494 throw DataException("Error - MeshAdapter::dump: copy Tags_keys to netCDF buffer failed: " + *newFileName);
495
496 // Tags_names_*
497 // This is an array of strings, it should be stored as an array but instead I have hacked in one attribute per string
498 // because the NetCDF manual doesn't tell how to do an array of strings
499 tag_map = mesh->TagMap;
500 if (tag_map) {
501 int i = 0;
502 while (tag_map) {
503 sprintf(name_temp, "Tags_name_%d", i);
504 if (!dataFile.add_att(name_temp, tag_map->name) )
505 throw DataException("Error - MeshAdapter::dump: appending Tags_names_ to NetCDF file failed: " + *newFileName);
506 tag_map=tag_map->next;
507 i++;
508 }
509 }
510
511 TMPMEMFREE(Tags_keys);
512
513 }
514
515
516 // NetCDF file is closed by destructor of NcFile object
517 #else
518 Finley_setError(IO_ERROR, "MeshAdapter::dump: not configured with NetCDF. Please contact your installation manager.");
519 #endif /* USE_NETCDF */
520 checkFinleyError();
521 }
522
523 string MeshAdapter::getDescription() const
524 {
525 return "FinleyMesh";
526 }
527
528 string MeshAdapter::functionSpaceTypeAsString(int functionSpaceType) const
529 {
530 FunctionSpaceNamesMapType::iterator loc;
531 loc=m_functionSpaceTypeNames.find(functionSpaceType);
532 if (loc==m_functionSpaceTypeNames.end()) {
533 return "Invalid function space type code.";
534 } else {
535 return loc->second;
536 }
537 }
538
539 bool MeshAdapter::isValidFunctionSpaceType(int functionSpaceType) const
540 {
541 FunctionSpaceNamesMapType::iterator loc;
542 loc=m_functionSpaceTypeNames.find(functionSpaceType);
543 return (loc!=m_functionSpaceTypeNames.end());
544 }
545
546 void MeshAdapter::setFunctionSpaceTypeNames()
547 {
548 m_functionSpaceTypeNames.insert
549 (FunctionSpaceNamesMapType::value_type(DegreesOfFreedom,"Finley_DegreesOfFreedom"));
550 m_functionSpaceTypeNames.insert
551 (FunctionSpaceNamesMapType::value_type(ReducedDegreesOfFreedom,"Finley_ReducedDegreesOfFreedom"));
552 m_functionSpaceTypeNames.insert
553 (FunctionSpaceNamesMapType::value_type(Nodes,"Finley_Nodes"));
554 m_functionSpaceTypeNames.insert
555 (FunctionSpaceNamesMapType::value_type(ReducedNodes,"Finley_Reduced_Nodes"));
556 m_functionSpaceTypeNames.insert
557 (FunctionSpaceNamesMapType::value_type(Elements,"Finley_Elements"));
558 m_functionSpaceTypeNames.insert
559 (FunctionSpaceNamesMapType::value_type(ReducedElements,"Finley_Reduced_Elements"));
560 m_functionSpaceTypeNames.insert
561 (FunctionSpaceNamesMapType::value_type(FaceElements,"Finley_Face_Elements"));
562 m_functionSpaceTypeNames.insert
563 (FunctionSpaceNamesMapType::value_type(ReducedFaceElements,"Finley_Reduced_Face_Elements"));
564 m_functionSpaceTypeNames.insert
565 (FunctionSpaceNamesMapType::value_type(Points,"Finley_Points"));
566 m_functionSpaceTypeNames.insert
567 (FunctionSpaceNamesMapType::value_type(ContactElementsZero,"Finley_Contact_Elements_0"));
568 m_functionSpaceTypeNames.insert
569 (FunctionSpaceNamesMapType::value_type(ReducedContactElementsZero,"Finley_Reduced_Contact_Elements_0"));
570 m_functionSpaceTypeNames.insert
571 (FunctionSpaceNamesMapType::value_type(ContactElementsOne,"Finley_Contact_Elements_1"));
572 m_functionSpaceTypeNames.insert
573 (FunctionSpaceNamesMapType::value_type(ReducedContactElementsOne,"Finley_Reduced_Contact_Elements_1"));
574 }
575
576 int MeshAdapter::getContinuousFunctionCode() const
577 {
578 return Nodes;
579 }
580 int MeshAdapter::getReducedContinuousFunctionCode() const
581 {
582 return ReducedNodes;
583 }
584
585 int MeshAdapter::getFunctionCode() const
586 {
587 return Elements;
588 }
589 int MeshAdapter::getReducedFunctionCode() const
590 {
591 return ReducedElements;
592 }
593
594 int MeshAdapter::getFunctionOnBoundaryCode() const
595 {
596 return FaceElements;
597 }
598 int MeshAdapter::getReducedFunctionOnBoundaryCode() const
599 {
600 return ReducedFaceElements;
601 }
602
603 int MeshAdapter::getFunctionOnContactZeroCode() const
604 {
605 return ContactElementsZero;
606 }
607 int MeshAdapter::getReducedFunctionOnContactZeroCode() const
608 {
609 return ReducedContactElementsZero;
610 }
611
612 int MeshAdapter::getFunctionOnContactOneCode() const
613 {
614 return ContactElementsOne;
615 }
616 int MeshAdapter::getReducedFunctionOnContactOneCode() const
617 {
618 return ReducedContactElementsOne;
619 }
620
621 int MeshAdapter::getSolutionCode() const
622 {
623 return DegreesOfFreedom;
624 }
625
626 int MeshAdapter::getReducedSolutionCode() const
627 {
628 return ReducedDegreesOfFreedom;
629 }
630
631 int MeshAdapter::getDiracDeltaFunctionCode() const
632 {
633 return Points;
634 }
635
636 //
637 // return the spatial dimension of the Mesh:
638 //
639 int MeshAdapter::getDim() const
640 {
641 int numDim=Finley_Mesh_getDim(m_finleyMesh.get());
642 checkFinleyError();
643 return numDim;
644 }
645
646 //
647 // return the number of data points per sample and the number of samples
648 // needed to represent data on a parts of the mesh.
649 //
650 pair<int,int> MeshAdapter::getDataShape(int functionSpaceCode) const
651 {
652 int numDataPointsPerSample=0;
653 int numSamples=0;
654 Finley_Mesh* mesh=m_finleyMesh.get();
655 switch (functionSpaceCode) {
656 case(Nodes):
657 numDataPointsPerSample=1;
658 numSamples=Finley_NodeFile_getNumNodes(mesh->Nodes);
659 break;
660 case(ReducedNodes):
661 numDataPointsPerSample=1;
662 numSamples=Finley_NodeFile_getNumReducedNodes(mesh->Nodes);
663 break;
664 case(Elements):
665 if (mesh->Elements!=NULL) {
666 numSamples=mesh->Elements->numElements;
667 numDataPointsPerSample=mesh->Elements->ReferenceElement->numQuadNodes;
668 }
669 break;
670 case(ReducedElements):
671 if (mesh->Elements!=NULL) {
672 numSamples=mesh->Elements->numElements;
673 numDataPointsPerSample=mesh->Elements->ReferenceElementReducedOrder->numQuadNodes;
674 }
675 break;
676 case(FaceElements):
677 if (mesh->FaceElements!=NULL) {
678 numDataPointsPerSample=mesh->FaceElements->ReferenceElement->numQuadNodes;
679 numSamples=mesh->FaceElements->numElements;
680 }
681 break;
682 case(ReducedFaceElements):
683 if (mesh->FaceElements!=NULL) {
684 numDataPointsPerSample=mesh->FaceElements->ReferenceElementReducedOrder->numQuadNodes;
685 numSamples=mesh->FaceElements->numElements;
686 }
687 break;
688 case(Points):
689 if (mesh->Points!=NULL) {
690 numDataPointsPerSample=1;
691 numSamples=mesh->Points->numElements;
692 }
693 break;
694 case(ContactElementsZero):
695 if (mesh->ContactElements!=NULL) {
696 numDataPointsPerSample=mesh->ContactElements->ReferenceElement->numQuadNodes;
697 numSamples=mesh->ContactElements->numElements;
698 }
699 break;
700 case(ReducedContactElementsZero):
701 if (mesh->ContactElements!=NULL) {
702 numDataPointsPerSample=mesh->ContactElements->ReferenceElementReducedOrder->numQuadNodes;
703 numSamples=mesh->ContactElements->numElements;
704 }
705 break;
706 case(ContactElementsOne):
707 if (mesh->ContactElements!=NULL) {
708 numDataPointsPerSample=mesh->ContactElements->ReferenceElement->numQuadNodes;
709 numSamples=mesh->ContactElements->numElements;
710 }
711 break;
712 case(ReducedContactElementsOne):
713 if (mesh->ContactElements!=NULL) {
714 numDataPointsPerSample=mesh->ContactElements->ReferenceElementReducedOrder->numQuadNodes;
715 numSamples=mesh->ContactElements->numElements;
716 }
717 break;
718 case(DegreesOfFreedom):
719 if (mesh->Nodes!=NULL) {
720 numDataPointsPerSample=1;
721 numSamples=Finley_NodeFile_getNumDegreesOfFreedom(mesh->Nodes);
722 }
723 break;
724 case(ReducedDegreesOfFreedom):
725 if (mesh->Nodes!=NULL) {
726 numDataPointsPerSample=1;
727 numSamples=Finley_NodeFile_getNumReducedDegreesOfFreedom(mesh->Nodes);
728 }
729 break;
730 default:
731 stringstream temp;
732 temp << "Error - Invalid function space type: " << functionSpaceCode << " for domain: " << getDescription();
733 throw FinleyAdapterException(temp.str());
734 break;
735 }
736 return pair<int,int>(numDataPointsPerSample,numSamples);
737 }
738
739 //
740 // adds linear PDE of second order into a given stiffness matrix and right hand side:
741 //
742 void MeshAdapter::addPDEToSystem(
743 SystemMatrixAdapter& mat, escript::Data& rhs,
744 const escript::Data& A, const escript::Data& B, const escript::Data& C,const escript::Data& D,const escript::Data& X,const escript::Data& Y,
745 const escript::Data& d, const escript::Data& y,
746 const escript::Data& d_contact,const escript::Data& y_contact) const
747 {
748 escriptDataC _rhs=rhs.getDataC();
749 escriptDataC _A =A.getDataC();
750 escriptDataC _B=B.getDataC();
751 escriptDataC _C=C.getDataC();
752 escriptDataC _D=D.getDataC();
753 escriptDataC _X=X.getDataC();
754 escriptDataC _Y=Y.getDataC();
755 escriptDataC _d=d.getDataC();
756 escriptDataC _y=y.getDataC();
757 escriptDataC _d_contact=d_contact.getDataC();
758 escriptDataC _y_contact=y_contact.getDataC();
759
760 Finley_Mesh* mesh=m_finleyMesh.get();
761
762 Finley_Assemble_PDE(mesh->Nodes,mesh->Elements,mat.getPaso_SystemMatrix(), &_rhs, &_A, &_B, &_C, &_D, &_X, &_Y );
763 checkFinleyError();
764
765 Finley_Assemble_PDE(mesh->Nodes,mesh->FaceElements, mat.getPaso_SystemMatrix(), &_rhs, 0, 0, 0, &_d, 0, &_y );
766 checkFinleyError();
767
768 Finley_Assemble_PDE(mesh->Nodes,mesh->ContactElements, mat.getPaso_SystemMatrix(), &_rhs , 0, 0, 0, &_d_contact, 0, &_y_contact );
769 checkFinleyError();
770 }
771
772 void MeshAdapter::addPDEToLumpedSystem(
773 escript::Data& mat,
774 const escript::Data& D,
775 const escript::Data& d) const
776 {
777 escriptDataC _mat=mat.getDataC();
778 escriptDataC _D=D.getDataC();
779 escriptDataC _d=d.getDataC();
780
781 Finley_Mesh* mesh=m_finleyMesh.get();
782
783 Finley_Assemble_LumpedSystem(mesh->Nodes,mesh->Elements,&_mat, &_D);
784 Finley_Assemble_LumpedSystem(mesh->Nodes,mesh->FaceElements,&_mat, &_d);
785
786 checkFinleyError();
787 }
788
789
790 //
791 // adds linear PDE of second order into the right hand side only
792 //
793 void MeshAdapter::addPDEToRHS( escript::Data& rhs, const escript::Data& X,const escript::Data& Y, const escript::Data& y, const escript::Data& y_contact) const
794 {
795 Finley_Mesh* mesh=m_finleyMesh.get();
796
797 escriptDataC _rhs=rhs.getDataC();
798 escriptDataC _X=X.getDataC();
799 escriptDataC _Y=Y.getDataC();
800 escriptDataC _y=y.getDataC();
801 escriptDataC _y_contact=y_contact.getDataC();
802
803 Finley_Assemble_PDE(mesh->Nodes,mesh->Elements, 0, &_rhs, 0, 0, 0, 0, &_X, &_Y );
804 checkFinleyError();
805
806 Finley_Assemble_PDE(mesh->Nodes,mesh->FaceElements, 0, &_rhs, 0, 0, 0, 0, 0, &_y );
807 checkFinleyError();
808
809 Finley_Assemble_PDE(mesh->Nodes,mesh->ContactElements, 0, &_rhs , 0, 0, 0, 0, 0, &_y_contact );
810 checkFinleyError();
811 }
812
813 //
814 // interpolates data between different function spaces:
815 //
816 void MeshAdapter::interpolateOnDomain(escript::Data& target,const escript::Data& in) const
817 {
818 const MeshAdapter& inDomain=dynamic_cast<const MeshAdapter&>(in.getFunctionSpace().getDomain());
819 const MeshAdapter& targetDomain=dynamic_cast<const MeshAdapter&>(target.getFunctionSpace().getDomain());
820 if (inDomain!=*this)
821 throw FinleyAdapterException("Error - Illegal domain of interpolant.");
822 if (targetDomain!=*this)
823 throw FinleyAdapterException("Error - Illegal domain of interpolation target.");
824
825 Finley_Mesh* mesh=m_finleyMesh.get();
826 escriptDataC _target=target.getDataC();
827 escriptDataC _in=in.getDataC();
828 switch(in.getFunctionSpace().getTypeCode()) {
829 case(Nodes):
830 switch(target.getFunctionSpace().getTypeCode()) {
831 case(Nodes):
832 case(ReducedNodes):
833 case(DegreesOfFreedom):
834 case(ReducedDegreesOfFreedom):
835 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
836 break;
837 case(Elements):
838 case(ReducedElements):
839 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in,&_target);
840 break;
841 case(FaceElements):
842 case(ReducedFaceElements):
843 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in,&_target);
844 break;
845 case(Points):
846 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in,&_target);
847 break;
848 case(ContactElementsZero):
849 case(ReducedContactElementsZero):
850 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
851 break;
852 case(ContactElementsOne):
853 case(ReducedContactElementsOne):
854 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
855 break;
856 default:
857 stringstream temp;
858 temp << "Error - Interpolation on Domain: Finley does not know anything about function space type " << target.getFunctionSpace().getTypeCode();
859 throw FinleyAdapterException(temp.str());
860 break;
861 }
862 break;
863 case(ReducedNodes):
864 switch(target.getFunctionSpace().getTypeCode()) {
865 case(Nodes):
866 case(ReducedNodes):
867 case(DegreesOfFreedom):
868 case(ReducedDegreesOfFreedom):
869 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
870 break;
871 case(Elements):
872 case(ReducedElements):
873 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in,&_target);
874 break;
875 case(FaceElements):
876 case(ReducedFaceElements):
877 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in,&_target);
878 break;
879 case(Points):
880 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in,&_target);
881 break;
882 case(ContactElementsZero):
883 case(ReducedContactElementsZero):
884 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
885 break;
886 case(ContactElementsOne):
887 case(ReducedContactElementsOne):
888 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
889 break;
890 default:
891 stringstream temp;
892 temp << "Error - Interpolation on Domain: Finley does not know anything about function space type " << target.getFunctionSpace().getTypeCode();
893 throw FinleyAdapterException(temp.str());
894 break;
895 }
896 break;
897 case(Elements):
898 if (target.getFunctionSpace().getTypeCode()==Elements) {
899 Finley_Assemble_CopyElementData(mesh->Elements,&_target,&_in);
900 } else if (target.getFunctionSpace().getTypeCode()==ReducedElements) {
901 Finley_Assemble_AverageElementData(mesh->Elements,&_target,&_in);
902 } else {
903 throw FinleyAdapterException("Error - No interpolation with data on elements possible.");
904 }
905 break;
906 case(ReducedElements):
907 if (target.getFunctionSpace().getTypeCode()==ReducedElements) {
908 Finley_Assemble_CopyElementData(mesh->Elements,&_target,&_in);
909 } else {
910 throw FinleyAdapterException("Error - No interpolation with data on elements with reduced integration order possible.");
911 }
912 break;
913 case(FaceElements):
914 if (target.getFunctionSpace().getTypeCode()==FaceElements) {
915 Finley_Assemble_CopyElementData(mesh->FaceElements,&_target,&_in);
916 } else if (target.getFunctionSpace().getTypeCode()==ReducedFaceElements) {
917 Finley_Assemble_AverageElementData(mesh->FaceElements,&_target,&_in);
918 } else {
919 throw FinleyAdapterException("Error - No interpolation with data on face elements possible.");
920 }
921 break;
922 case(ReducedFaceElements):
923 if (target.getFunctionSpace().getTypeCode()==FaceElements) {
924 Finley_Assemble_CopyElementData(mesh->FaceElements,&_target,&_in);
925 } else {
926 throw FinleyAdapterException("Error - No interpolation with data on face elements with reduced integration order possible.");
927 }
928 break;
929 case(Points):
930 if (target.getFunctionSpace().getTypeCode()==Points) {
931 Finley_Assemble_CopyElementData(mesh->Points,&_target,&_in);
932 } else {
933 throw FinleyAdapterException("Error - No interpolation with data on points possible.");
934 }
935 break;
936 case(ContactElementsZero):
937 case(ContactElementsOne):
938 if (target.getFunctionSpace().getTypeCode()==ContactElementsZero || target.getFunctionSpace().getTypeCode()==ContactElementsOne) {
939 Finley_Assemble_CopyElementData(mesh->ContactElements,&_target,&_in);
940 } else if (target.getFunctionSpace().getTypeCode()==ReducedContactElementsZero || target.getFunctionSpace().getTypeCode()==ReducedContactElementsOne) {
941 Finley_Assemble_AverageElementData(mesh->ContactElements,&_target,&_in);
942 } else {
943 throw FinleyAdapterException("Error - No interpolation with data on contact elements possible.");
944 break;
945 }
946 break;
947 case(ReducedContactElementsZero):
948 case(ReducedContactElementsOne):
949 if (target.getFunctionSpace().getTypeCode()==ReducedContactElementsZero || target.getFunctionSpace().getTypeCode()==ReducedContactElementsOne) {
950 Finley_Assemble_CopyElementData(mesh->ContactElements,&_target,&_in);
951 } else {
952 throw FinleyAdapterException("Error - No interpolation with data on contact elements with reduced integration order possible.");
953 break;
954 }
955 break;
956 case(DegreesOfFreedom):
957 switch(target.getFunctionSpace().getTypeCode()) {
958 case(ReducedDegreesOfFreedom):
959 case(DegreesOfFreedom):
960 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
961 break;
962
963 case(Nodes):
964 case(ReducedNodes):
965 if (getMPISize()>1) {
966 escript::Data temp=escript::Data(in);
967 temp.expand();
968 escriptDataC _in2 = temp.getDataC();
969 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in2);
970 } else {
971 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
972 }
973 break;
974 case(Elements):
975 case(ReducedElements):
976 if (getMPISize()>1) {
977 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
978 escriptDataC _in2 = temp.getDataC();
979 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in2,&_target);
980 } else {
981 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in,&_target);
982 }
983 break;
984 case(FaceElements):
985 case(ReducedFaceElements):
986 if (getMPISize()>1) {
987 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
988 escriptDataC _in2 = temp.getDataC();
989 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in2,&_target);
990
991 } else {
992 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in,&_target);
993 }
994 break;
995 case(Points):
996 if (getMPISize()>1) {
997 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
998 escriptDataC _in2 = temp.getDataC();
999 } else {
1000 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in,&_target);
1001 }
1002 break;
1003 case(ContactElementsZero):
1004 case(ContactElementsOne):
1005 case(ReducedContactElementsZero):
1006 case(ReducedContactElementsOne):
1007 if (getMPISize()>1) {
1008 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
1009 escriptDataC _in2 = temp.getDataC();
1010 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in2,&_target);
1011 } else {
1012 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
1013 }
1014 break;
1015 default:
1016 stringstream temp;
1017 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << target.getFunctionSpace().getTypeCode();
1018 throw FinleyAdapterException(temp.str());
1019 break;
1020 }
1021 break;
1022 case(ReducedDegreesOfFreedom):
1023 switch(target.getFunctionSpace().getTypeCode()) {
1024 case(Nodes):
1025 throw FinleyAdapterException("Error - Finley does not support interpolation from reduced degrees of freedom to mesh nodes.");
1026 break;
1027 case(ReducedNodes):
1028 if (getMPISize()>1) {
1029 escript::Data temp=escript::Data(in);
1030 temp.expand();
1031 escriptDataC _in2 = temp.getDataC();
1032 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in2);
1033 } else {
1034 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
1035 }
1036 break;
1037 case(DegreesOfFreedom):
1038 throw FinleyAdapterException("Error - Finley does not support interpolation from reduced degrees of freedom to degrees of freedom");
1039 break;
1040 case(ReducedDegreesOfFreedom):
1041 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
1042 break;
1043 case(Elements):
1044 case(ReducedElements):
1045 if (getMPISize()>1) {
1046 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
1047 escriptDataC _in2 = temp.getDataC();
1048 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in2,&_target);
1049 } else {
1050 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in,&_target);
1051 }
1052 break;
1053 case(FaceElements):
1054 case(ReducedFaceElements):
1055 if (getMPISize()>1) {
1056 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
1057 escriptDataC _in2 = temp.getDataC();
1058 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in2,&_target);
1059 } else {
1060 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in,&_target);
1061 }
1062 break;
1063 case(Points):
1064 if (getMPISize()>1) {
1065 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
1066 escriptDataC _in2 = temp.getDataC();
1067 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in2,&_target);
1068 } else {
1069 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in,&_target);
1070 }
1071 break;
1072 case(ContactElementsZero):
1073 case(ContactElementsOne):
1074 case(ReducedContactElementsZero):
1075 case(ReducedContactElementsOne):
1076 if (getMPISize()>1) {
1077 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
1078 escriptDataC _in2 = temp.getDataC();
1079 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in2,&_target);
1080 } else {
1081 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
1082 }
1083 break;
1084 default:
1085 stringstream temp;
1086 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << target.getFunctionSpace().getTypeCode();
1087 throw FinleyAdapterException(temp.str());
1088 break;
1089 }
1090 break;
1091 default:
1092 stringstream temp;
1093 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type %d" << in.getFunctionSpace().getTypeCode();
1094 throw FinleyAdapterException(temp.str());
1095 break;
1096 }
1097 checkFinleyError();
1098 }
1099
1100 //
1101 // copies the locations of sample points into x:
1102 //
1103 void MeshAdapter::setToX(escript::Data& arg) const
1104 {
1105 const MeshAdapter& argDomain=dynamic_cast<const MeshAdapter&>(arg.getFunctionSpace().getDomain());
1106 if (argDomain!=*this)
1107 throw FinleyAdapterException("Error - Illegal domain of data point locations");
1108 Finley_Mesh* mesh=m_finleyMesh.get();
1109 // in case of values node coordinates we can do the job directly:
1110 if (arg.getFunctionSpace().getTypeCode()==Nodes) {
1111 escriptDataC _arg=arg.getDataC();
1112 Finley_Assemble_NodeCoordinates(mesh->Nodes,&_arg);
1113 } else {
1114 escript::Data tmp_data=Vector(0.0,continuousFunction(asAbstractContinuousDomain()),true);
1115 escriptDataC _tmp_data=tmp_data.getDataC();
1116 Finley_Assemble_NodeCoordinates(mesh->Nodes,&_tmp_data);
1117 // this is then interpolated onto arg:
1118 interpolateOnDomain(arg,tmp_data);
1119 }
1120 checkFinleyError();
1121 }
1122
1123 //
1124 // return the normal vectors at the location of data points as a Data object:
1125 //
1126 void MeshAdapter::setToNormal(escript::Data& normal) const
1127 {
1128 const MeshAdapter& normalDomain=dynamic_cast<const MeshAdapter&>(normal.getFunctionSpace().getDomain());
1129 if (normalDomain!=*this)
1130 throw FinleyAdapterException("Error - Illegal domain of normal locations");
1131 Finley_Mesh* mesh=m_finleyMesh.get();
1132 escriptDataC _normal=normal.getDataC();
1133 switch(normal.getFunctionSpace().getTypeCode()) {
1134 case(Nodes):
1135 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for nodes");
1136 break;
1137 case(ReducedNodes):
1138 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for reduced nodes");
1139 break;
1140 case(Elements):
1141 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for elements");
1142 break;
1143 case(ReducedElements):
1144 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for elements with reduced integration order");
1145 break;
1146 case (FaceElements):
1147 Finley_Assemble_setNormal(mesh->Nodes,mesh->FaceElements,&_normal);
1148 break;
1149 case (ReducedFaceElements):
1150 Finley_Assemble_setNormal(mesh->Nodes,mesh->FaceElements,&_normal);
1151 break;
1152 case(Points):
1153 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for point elements");
1154 break;
1155 case (ContactElementsOne):
1156 case (ContactElementsZero):
1157 Finley_Assemble_setNormal(mesh->Nodes,mesh->ContactElements,&_normal);
1158 break;
1159 case (ReducedContactElementsOne):
1160 case (ReducedContactElementsZero):
1161 Finley_Assemble_setNormal(mesh->Nodes,mesh->ContactElements,&_normal);
1162 break;
1163 case(DegreesOfFreedom):
1164 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for degrees of freedom.");
1165 break;
1166 case(ReducedDegreesOfFreedom):
1167 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for reduced degrees of freedom.");
1168 break;
1169 default:
1170 stringstream temp;
1171 temp << "Error - Normal Vectors: Finley does not know anything about function space type " << normal.getFunctionSpace().getTypeCode();
1172 throw FinleyAdapterException(temp.str());
1173 break;
1174 }
1175 checkFinleyError();
1176 }
1177
1178 //
1179 // interpolates data to other domain:
1180 //
1181 void MeshAdapter::interpolateACross(escript::Data& target,const escript::Data& source) const
1182 {
1183 const MeshAdapter& targetDomain=dynamic_cast<const MeshAdapter&>(target.getFunctionSpace().getDomain());
1184 if (targetDomain!=*this)
1185 throw FinleyAdapterException("Error - Illegal domain of interpolation target");
1186
1187 throw FinleyAdapterException("Error - Finley does not allow interpolation across domains yet.");
1188 }
1189
1190 //
1191 // calculates the integral of a function defined of arg:
1192 //
1193 void MeshAdapter::setToIntegrals(std::vector<double>& integrals,const escript::Data& arg) const
1194 {
1195 const MeshAdapter& argDomain=dynamic_cast<const MeshAdapter&>(arg.getFunctionSpace().getDomain());
1196 if (argDomain!=*this)
1197 throw FinleyAdapterException("Error - Illegal domain of integration kernel");
1198
1199 double blocktimer_start = blocktimer_time();
1200 Finley_Mesh* mesh=m_finleyMesh.get();
1201 escriptDataC _temp;
1202 escript::Data temp;
1203 escriptDataC _arg=arg.getDataC();
1204 switch(arg.getFunctionSpace().getTypeCode()) {
1205 case(Nodes):
1206 temp=escript::Data( arg, function(asAbstractContinuousDomain()) );
1207 _temp=temp.getDataC();
1208 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_temp,&integrals[0]);
1209 break;
1210 case(ReducedNodes):
1211 temp=escript::Data( arg, function(asAbstractContinuousDomain()) );
1212 _temp=temp.getDataC();
1213 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_temp,&integrals[0]);
1214 break;
1215 case(Elements):
1216 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_arg,&integrals[0]);
1217 break;
1218 case(ReducedElements):
1219 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_arg,&integrals[0]);
1220 break;
1221 case(FaceElements):
1222 Finley_Assemble_integrate(mesh->Nodes,mesh->FaceElements,&_arg,&integrals[0]);
1223 break;
1224 case(ReducedFaceElements):
1225 Finley_Assemble_integrate(mesh->Nodes,mesh->FaceElements,&_arg,&integrals[0]);
1226 break;
1227 case(Points):
1228 throw FinleyAdapterException("Error - Integral of data on points is not supported.");
1229 break;
1230 case(ContactElementsZero):
1231 Finley_Assemble_integrate(mesh->Nodes,mesh->ContactElements,&_arg,&integrals[0]);
1232 break;
1233 case(ReducedContactElementsZero):
1234 Finley_Assemble_integrate(mesh->Nodes,mesh->ContactElements,&_arg,&integrals[0]);
1235 break;
1236 case(ContactElementsOne):
1237 Finley_Assemble_integrate(mesh->Nodes,mesh->ContactElements,&_arg,&integrals[0]);
1238 break;
1239 case(ReducedContactElementsOne):
1240 Finley_Assemble_integrate(mesh->Nodes,mesh->ContactElements,&_arg,&integrals[0]);
1241 break;
1242 case(DegreesOfFreedom):
1243 temp=escript::Data( arg, function(asAbstractContinuousDomain()) );
1244 _temp=temp.getDataC();
1245 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_temp,&integrals[0]);
1246 break;
1247 case(ReducedDegreesOfFreedom):
1248 temp=escript::Data( arg, function(asAbstractContinuousDomain()) );
1249 _temp=temp.getDataC();
1250 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_temp,&integrals[0]);
1251 break;
1252 default:
1253 stringstream temp;
1254 temp << "Error - Integrals: Finley does not know anything about function space type " << arg.getFunctionSpace().getTypeCode();
1255 throw FinleyAdapterException(temp.str());
1256 break;
1257 }
1258 checkFinleyError();
1259 blocktimer_increment("integrate()", blocktimer_start);
1260 }
1261
1262 //
1263 // calculates the gradient of arg:
1264 //
1265 void MeshAdapter::setToGradient(escript::Data& grad,const escript::Data& arg) const
1266 {
1267 const MeshAdapter& argDomain=dynamic_cast<const MeshAdapter&>(arg.getFunctionSpace().getDomain());
1268 if (argDomain!=*this)
1269 throw FinleyAdapterException("Error - Illegal domain of gradient argument");
1270 const MeshAdapter& gradDomain=dynamic_cast<const MeshAdapter&>(grad.getFunctionSpace().getDomain());
1271 if (gradDomain!=*this)
1272 throw FinleyAdapterException("Error - Illegal domain of gradient");
1273
1274 Finley_Mesh* mesh=m_finleyMesh.get();
1275 escriptDataC _grad=grad.getDataC();
1276 escriptDataC nodeDataC;
1277 escript::Data temp;
1278 if (getMPISize()>1) {
1279 if( arg.getFunctionSpace().getTypeCode() == DegreesOfFreedom ) {
1280 temp=escript::Data( arg, continuousFunction(asAbstractContinuousDomain()) );
1281 nodeDataC = temp.getDataC();
1282 } else if( arg.getFunctionSpace().getTypeCode() == ReducedDegreesOfFreedom ) {
1283 temp=escript::Data( arg, reducedContinuousFunction(asAbstractContinuousDomain()) );
1284 nodeDataC = temp.getDataC();
1285 } else {
1286 nodeDataC = arg.getDataC();
1287 }
1288 } else {
1289 nodeDataC = arg.getDataC();
1290 }
1291 switch(grad.getFunctionSpace().getTypeCode()) {
1292 case(Nodes):
1293 throw FinleyAdapterException("Error - Gradient at nodes is not supported.");
1294 break;
1295 case(ReducedNodes):
1296 throw FinleyAdapterException("Error - Gradient at reduced nodes is not supported.");
1297 break;
1298 case(Elements):
1299 Finley_Assemble_gradient(mesh->Nodes,mesh->Elements,&_grad,&nodeDataC);
1300 break;
1301 case(ReducedElements):
1302 Finley_Assemble_gradient(mesh->Nodes,mesh->Elements,&_grad,&nodeDataC);
1303 break;
1304 case(FaceElements):
1305 Finley_Assemble_gradient(mesh->Nodes,mesh->FaceElements,&_grad,&nodeDataC);
1306 break;
1307 case(ReducedFaceElements):
1308 Finley_Assemble_gradient(mesh->Nodes,mesh->FaceElements,&_grad,&nodeDataC);
1309 break;
1310 case(Points):
1311 throw FinleyAdapterException("Error - Gradient at points is not supported.");
1312 break;
1313 case(ContactElementsZero):
1314 Finley_Assemble_gradient(mesh->Nodes,mesh->ContactElements,&_grad,&nodeDataC);
1315 break;
1316 case(ReducedContactElementsZero):
1317 Finley_Assemble_gradient(mesh->Nodes,mesh->ContactElements,&_grad,&nodeDataC);
1318 break;
1319 case(ContactElementsOne):
1320 Finley_Assemble_gradient(mesh->Nodes,mesh->ContactElements,&_grad,&nodeDataC);
1321 break;
1322 case(ReducedContactElementsOne):
1323 Finley_Assemble_gradient(mesh->Nodes,mesh->ContactElements,&_grad,&nodeDataC);
1324 break;
1325 case(DegreesOfFreedom):
1326 throw FinleyAdapterException("Error - Gradient at degrees of freedom is not supported.");
1327 break;
1328 case(ReducedDegreesOfFreedom):
1329 throw FinleyAdapterException("Error - Gradient at reduced degrees of freedom is not supported.");
1330 break;
1331 default:
1332 stringstream temp;
1333 temp << "Error - Gradient: Finley does not know anything about function space type " << arg.getFunctionSpace().getTypeCode();
1334 throw FinleyAdapterException(temp.str());
1335 break;
1336 }
1337 checkFinleyError();
1338 }
1339
1340 //
1341 // returns the size of elements:
1342 //
1343 void MeshAdapter::setToSize(escript::Data& size) const
1344 {
1345 Finley_Mesh* mesh=m_finleyMesh.get();
1346 escriptDataC tmp=size.getDataC();
1347 switch(size.getFunctionSpace().getTypeCode()) {
1348 case(Nodes):
1349 throw FinleyAdapterException("Error - Size of nodes is not supported.");
1350 break;
1351 case(ReducedNodes):
1352 throw FinleyAdapterException("Error - Size of reduced nodes is not supported.");
1353 break;
1354 case(Elements):
1355 Finley_Assemble_getSize(mesh->Nodes,mesh->Elements,&tmp);
1356 break;
1357 case(ReducedElements):
1358 Finley_Assemble_getSize(mesh->Nodes,mesh->Elements,&tmp);
1359 break;
1360 case(FaceElements):
1361 Finley_Assemble_getSize(mesh->Nodes,mesh->FaceElements,&tmp);
1362 break;
1363 case(ReducedFaceElements):
1364 Finley_Assemble_getSize(mesh->Nodes,mesh->FaceElements,&tmp);
1365 break;
1366 case(Points):
1367 throw FinleyAdapterException("Error - Size of point elements is not supported.");
1368 break;
1369 case(ContactElementsZero):
1370 case(ContactElementsOne):
1371 Finley_Assemble_getSize(mesh->Nodes,mesh->ContactElements,&tmp);
1372 break;
1373 case(ReducedContactElementsZero):
1374 case(ReducedContactElementsOne):
1375 Finley_Assemble_getSize(mesh->Nodes,mesh->ContactElements,&tmp);
1376 break;
1377 case(DegreesOfFreedom):
1378 throw FinleyAdapterException("Error - Size of degrees of freedom is not supported.");
1379 break;
1380 case(ReducedDegreesOfFreedom):
1381 throw FinleyAdapterException("Error - Size of reduced degrees of freedom is not supported.");
1382 break;
1383 default:
1384 stringstream temp;
1385 temp << "Error - Element size: Finley does not know anything about function space type " << size.getFunctionSpace().getTypeCode();
1386 throw FinleyAdapterException(temp.str());
1387 break;
1388 }
1389 checkFinleyError();
1390 }
1391
1392 // sets the location of nodes:
1393 void MeshAdapter::setNewX(const escript::Data& new_x)
1394 {
1395 Finley_Mesh* mesh=m_finleyMesh.get();
1396 escriptDataC tmp;
1397 const MeshAdapter& newDomain=dynamic_cast<const MeshAdapter&>(new_x.getFunctionSpace().getDomain());
1398 if (newDomain!=*this)
1399 throw FinleyAdapterException("Error - Illegal domain of new point locations");
1400 tmp = new_x.getDataC();
1401 Finley_Mesh_setCoordinates(mesh,&tmp);
1402 checkFinleyError();
1403 }
1404
1405 // saves a data array in openDX format:
1406 void MeshAdapter::saveDX(const std::string& filename,const boost::python::dict& arg) const
1407 {
1408 int MAX_namelength=256;
1409 const int num_data=boost::python::extract<int>(arg.attr("__len__")());
1410 /* win32 refactor */
1411 char* *names = (num_data>0) ? TMPMEMALLOC(num_data,char*) : (char**)NULL;
1412 for(int i=0;i<num_data;i++)
1413 {
1414 names[i] = (MAX_namelength>0) ? TMPMEMALLOC(MAX_namelength,char) : (char*)NULL;
1415 }
1416
1417 char* *c_names = (num_data>0) ? TMPMEMALLOC(num_data,char*) : (char**)NULL;
1418 escriptDataC *data = (num_data>0) ? TMPMEMALLOC(num_data,escriptDataC) : (escriptDataC*)NULL;
1419 escriptDataC* *ptr_data = (num_data>0) ? TMPMEMALLOC(num_data,escriptDataC*) : (escriptDataC**)NULL;
1420
1421 boost::python::list keys=arg.keys();
1422 for (int i=0;i<num_data;++i) {
1423 std::string n=boost::python::extract<std::string>(keys[i]);
1424 escript::Data& d=boost::python::extract<escript::Data&>(arg[keys[i]]);
1425 if (dynamic_cast<const MeshAdapter&>(d.getFunctionSpace().getDomain()) !=*this)
1426 throw FinleyAdapterException("Error in saveVTK: Data must be defined on same Domain");
1427 data[i]=d.getDataC();
1428 ptr_data[i]=&(data[i]);
1429 c_names[i]=&(names[i][0]);
1430 if (n.length()>MAX_namelength-1) {
1431 strncpy(c_names[i],n.c_str(),MAX_namelength-1);
1432 c_names[i][MAX_namelength-1]='\0';
1433 } else {
1434 strcpy(c_names[i],n.c_str());
1435 }
1436 }
1437 Finley_Mesh_saveDX(filename.c_str(),m_finleyMesh.get(),num_data,c_names,ptr_data);
1438 checkFinleyError();
1439
1440 /* win32 refactor */
1441 TMPMEMFREE(c_names);
1442 TMPMEMFREE(data);
1443 TMPMEMFREE(ptr_data);
1444 for(int i=0;i<num_data;i++)
1445 {
1446 TMPMEMFREE(names[i]);
1447 }
1448 TMPMEMFREE(names);
1449
1450 return;
1451 }
1452
1453 // saves a data array in openVTK format:
1454 void MeshAdapter::saveVTK(const std::string& filename,const boost::python::dict& arg) const
1455 {
1456 int MAX_namelength=256;
1457 const int num_data=boost::python::extract<int>(arg.attr("__len__")());
1458 /* win32 refactor */
1459 char* *names = (num_data>0) ? TMPMEMALLOC(num_data,char*) : (char**)NULL;
1460 for(int i=0;i<num_data;i++)
1461 {
1462 names[i] = (MAX_namelength>0) ? TMPMEMALLOC(MAX_namelength,char) : (char*)NULL;
1463 }
1464
1465 char* *c_names = (num_data>0) ? TMPMEMALLOC(num_data,char*) : (char**)NULL;
1466 escriptDataC *data = (num_data>0) ? TMPMEMALLOC(num_data,escriptDataC) : (escriptDataC*)NULL;
1467 escriptDataC* *ptr_data = (num_data>0) ? TMPMEMALLOC(num_data,escriptDataC*) : (escriptDataC**)NULL;
1468
1469 boost::python::list keys=arg.keys();
1470 for (int i=0;i<num_data;++i) {
1471 std::string n=boost::python::extract<std::string>(keys[i]);
1472 escript::Data& d=boost::python::extract<escript::Data&>(arg[keys[i]]);
1473 if (dynamic_cast<const MeshAdapter&>(d.getFunctionSpace().getDomain()) !=*this)
1474 throw FinleyAdapterException("Error in saveVTK: Data must be defined on same Domain");
1475 data[i]=d.getDataC();
1476 ptr_data[i]=&(data[i]);
1477 c_names[i]=&(names[i][0]);
1478 if (n.length()>MAX_namelength-1) {
1479 strncpy(c_names[i],n.c_str(),MAX_namelength-1);
1480 c_names[i][MAX_namelength-1]='\0';
1481 } else {
1482 strcpy(c_names[i],n.c_str());
1483 }
1484 }
1485 Finley_Mesh_saveVTK(filename.c_str(),m_finleyMesh.get(),num_data,c_names,ptr_data);
1486
1487 checkFinleyError();
1488 /* win32 refactor */
1489 TMPMEMFREE(c_names);
1490 TMPMEMFREE(data);
1491 TMPMEMFREE(ptr_data);
1492 for(int i=0;i<num_data;i++)
1493 {
1494 TMPMEMFREE(names[i]);
1495 }
1496 TMPMEMFREE(names);
1497
1498 return;
1499 }
1500
1501
1502 // creates a SystemMatrixAdapter stiffness matrix an initializes it with zeros:
1503 SystemMatrixAdapter MeshAdapter::newSystemMatrix(
1504 const int row_blocksize,
1505 const escript::FunctionSpace& row_functionspace,
1506 const int column_blocksize,
1507 const escript::FunctionSpace& column_functionspace,
1508 const int type) const
1509 {
1510 int reduceRowOrder=0;
1511 int reduceColOrder=0;
1512 // is the domain right?
1513 const MeshAdapter& row_domain=dynamic_cast<const MeshAdapter&>(row_functionspace.getDomain());
1514 if (row_domain!=*this)
1515 throw FinleyAdapterException("Error - domain of row function space does not match the domain of matrix generator.");
1516 const MeshAdapter& col_domain=dynamic_cast<const MeshAdapter&>(column_functionspace.getDomain());
1517 if (col_domain!=*this)
1518 throw FinleyAdapterException("Error - domain of columnn function space does not match the domain of matrix generator.");
1519 // is the function space type right
1520 if (row_functionspace.getTypeCode()==DegreesOfFreedom) {
1521 reduceRowOrder=0;
1522 } else if (row_functionspace.getTypeCode()==ReducedDegreesOfFreedom) {
1523 reduceRowOrder=1;
1524 } else {
1525 throw FinleyAdapterException("Error - illegal function space type for system matrix rows.");
1526 }
1527 if (column_functionspace.getTypeCode()==DegreesOfFreedom) {
1528 reduceColOrder=0;
1529 } else if (column_functionspace.getTypeCode()==ReducedDegreesOfFreedom) {
1530 reduceColOrder=1;
1531 } else {
1532 throw FinleyAdapterException("Error - illegal function space type for system matrix columns.");
1533 }
1534 // generate matrix:
1535
1536 Paso_SystemMatrixPattern* fsystemMatrixPattern=Finley_getPattern(getFinley_Mesh(),reduceRowOrder,reduceColOrder);
1537 checkFinleyError();
1538 Paso_SystemMatrix* fsystemMatrix;
1539 int trilinos = 0;
1540 if (trilinos) {
1541 #ifdef TRILINOS
1542 /* Allocation Epetra_VrbMatrix here */
1543 #endif
1544 }
1545 else {
1546 fsystemMatrix=Paso_SystemMatrix_alloc(type,fsystemMatrixPattern,row_blocksize,column_blocksize);
1547 }
1548 checkPasoError();
1549 Paso_SystemMatrixPattern_free(fsystemMatrixPattern);
1550 return SystemMatrixAdapter(fsystemMatrix,row_blocksize,row_functionspace,column_blocksize,column_functionspace);
1551 }
1552
1553 //
1554 // vtkObject MeshAdapter::createVtkObject() const
1555 // TODO:
1556 //
1557
1558 //
1559 // returns true if data at the atom_type is considered as being cell centered:
1560 bool MeshAdapter::isCellOriented(int functionSpaceCode) const
1561 {
1562 switch(functionSpaceCode) {
1563 case(Nodes):
1564 case(DegreesOfFreedom):
1565 case(ReducedDegreesOfFreedom):
1566 return false;
1567 break;
1568 case(Elements):
1569 case(FaceElements):
1570 case(Points):
1571 case(ContactElementsZero):
1572 case(ContactElementsOne):
1573 case(ReducedElements):
1574 case(ReducedFaceElements):
1575 case(ReducedContactElementsZero):
1576 case(ReducedContactElementsOne):
1577 return true;
1578 break;
1579 default:
1580 stringstream temp;
1581 temp << "Error - Cell: Finley does not know anything about function space type " << functionSpaceCode;
1582 throw FinleyAdapterException(temp.str());
1583 break;
1584 }
1585 checkFinleyError();
1586 return false;
1587 }
1588
1589 bool MeshAdapter::probeInterpolationOnDomain(int functionSpaceType_source,int functionSpaceType_target) const
1590 {
1591 switch(functionSpaceType_source) {
1592 case(Nodes):
1593 switch(functionSpaceType_target) {
1594 case(Nodes):
1595 case(ReducedNodes):
1596 case(ReducedDegreesOfFreedom):
1597 case(DegreesOfFreedom):
1598 case(Elements):
1599 case(ReducedElements):
1600 case(FaceElements):
1601 case(ReducedFaceElements):
1602 case(Points):
1603 case(ContactElementsZero):
1604 case(ReducedContactElementsZero):
1605 case(ContactElementsOne):
1606 case(ReducedContactElementsOne):
1607 return true;
1608 default:
1609 stringstream temp;
1610 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_target;
1611 throw FinleyAdapterException(temp.str());
1612 }
1613 break;
1614 case(ReducedNodes):
1615 switch(functionSpaceType_target) {
1616 case(ReducedNodes):
1617 case(ReducedDegreesOfFreedom):
1618 case(Elements):
1619 case(ReducedElements):
1620 case(FaceElements):
1621 case(ReducedFaceElements):
1622 case(Points):
1623 case(ContactElementsZero):
1624 case(ReducedContactElementsZero):
1625 case(ContactElementsOne):
1626 case(ReducedContactElementsOne):
1627 return true;
1628 case(Nodes):
1629 case(DegreesOfFreedom):
1630 return false;
1631 default:
1632 stringstream temp;
1633 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_target;
1634 throw FinleyAdapterException(temp.str());
1635 }
1636 break;
1637 case(Elements):
1638 if (functionSpaceType_target==Elements) {
1639 return true;
1640 } else if (functionSpaceType_target==ReducedElements) {
1641 return true;
1642 } else {
1643 return false;
1644 }
1645 case(ReducedElements):
1646 if (functionSpaceType_target==ReducedElements) {
1647 return true;
1648 } else {
1649 return false;
1650 }
1651 case(FaceElements):
1652 if (functionSpaceType_target==FaceElements) {
1653 return true;
1654 } else if (functionSpaceType_target==ReducedFaceElements) {
1655 return true;
1656 } else {
1657 return false;
1658 }
1659 case(ReducedFaceElements):
1660 if (functionSpaceType_target==ReducedFaceElements) {
1661 return true;
1662 } else {
1663 return false;
1664 }
1665 case(Points):
1666 if (functionSpaceType_target==Points) {
1667 return true;
1668 } else {
1669 return false;
1670 }
1671 case(ContactElementsZero):
1672 case(ContactElementsOne):
1673 if (functionSpaceType_target==ContactElementsZero || functionSpaceType_target==ContactElementsOne) {
1674 return true;
1675 } else if (functionSpaceType_target==ReducedContactElementsZero || functionSpaceType_target==ReducedContactElementsOne) {
1676 return true;
1677 } else {
1678 return false;
1679 }
1680 case(ReducedContactElementsZero):
1681 case(ReducedContactElementsOne):
1682 if (functionSpaceType_target==ReducedContactElementsZero || functionSpaceType_target==ReducedContactElementsOne) {
1683 return true;
1684 } else {
1685 return false;
1686 }
1687 case(DegreesOfFreedom):
1688 switch(functionSpaceType_target) {
1689 case(ReducedDegreesOfFreedom):
1690 case(DegreesOfFreedom):
1691 case(Nodes):
1692 case(ReducedNodes):
1693 case(Elements):
1694 case(ReducedElements):
1695 case(Points):
1696 case(FaceElements):
1697 case(ReducedFaceElements):
1698 case(ContactElementsZero):
1699 case(ReducedContactElementsZero):
1700 case(ContactElementsOne):
1701 case(ReducedContactElementsOne):
1702 return true;
1703 default:
1704 stringstream temp;
1705 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_target;
1706 throw FinleyAdapterException(temp.str());
1707 }
1708 break;
1709 case(ReducedDegreesOfFreedom):
1710 switch(functionSpaceType_target) {
1711 case(ReducedDegreesOfFreedom):
1712 case(ReducedNodes):
1713 case(Elements):
1714 case(ReducedElements):
1715 case(FaceElements):
1716 case(ReducedFaceElements):
1717 case(Points):
1718 case(ContactElementsZero):
1719 case(ReducedContactElementsZero):
1720 case(ContactElementsOne):
1721 case(ReducedContactElementsOne):
1722 return true;
1723 case(Nodes):
1724 case(DegreesOfFreedom):
1725 return false;
1726 default:
1727 stringstream temp;
1728 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_target;
1729 throw FinleyAdapterException(temp.str());
1730 }
1731 break;
1732 default:
1733 stringstream temp;
1734 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_source;
1735 throw FinleyAdapterException(temp.str());
1736 break;
1737 }
1738 checkFinleyError();
1739 return false;
1740 }
1741
1742 bool MeshAdapter::probeInterpolationACross(int functionSpaceType_source,const AbstractDomain& targetDomain, int functionSpaceType_target) const
1743 {
1744 return false;
1745 }
1746
1747 bool MeshAdapter::operator==(const AbstractDomain& other) const
1748 {
1749 const MeshAdapter* temp=dynamic_cast<const MeshAdapter*>(&other);
1750 if (temp!=0) {
1751 return (m_finleyMesh==temp->m_finleyMesh);
1752 } else {
1753 return false;
1754 }
1755 }
1756
1757 bool MeshAdapter::operator!=(const AbstractDomain& other) const
1758 {
1759 return !(operator==(other));
1760 }
1761
1762 int MeshAdapter::getSystemMatrixTypeId(const int solver, const int package, const bool symmetry) const
1763 {
1764 int out=Paso_SystemMatrix_getSystemMatrixTypeId(SystemMatrixAdapter::mapOptionToPaso(solver),SystemMatrixAdapter::mapOptionToPaso(package),symmetry?1:0);
1765 checkPasoError();
1766 return out;
1767 }
1768
1769 escript::Data MeshAdapter::getX() const
1770 {
1771 return continuousFunction(asAbstractContinuousDomain()).getX();
1772 }
1773
1774 escript::Data MeshAdapter::getNormal() const
1775 {
1776 return functionOnBoundary(asAbstractContinuousDomain()).getNormal();
1777 }
1778
1779 escript::Data MeshAdapter::getSize() const
1780 {
1781 return function(asAbstractContinuousDomain()).getSize();
1782 }
1783
1784 int* MeshAdapter::borrowSampleReferenceIDs(int functionSpaceType) const
1785 {
1786 int *out=0,i;
1787 Finley_Mesh* mesh=m_finleyMesh.get();
1788 switch (functionSpaceType) {
1789 case(Nodes):
1790 out=mesh->Nodes->Id;
1791 break;
1792 case(ReducedNodes):
1793 out=mesh->Nodes->reducedNodesId;
1794 break;
1795 case(Elements):
1796 out=mesh->Elements->Id;
1797 break;
1798 case(ReducedElements):
1799 out=mesh->Elements->Id;
1800 break;
1801 case(FaceElements):
1802 out=mesh->FaceElements->Id;
1803 break;
1804 case(ReducedFaceElements):
1805 out=mesh->FaceElements->Id;
1806 break;
1807 case(Points):
1808 out=mesh->Points->Id;
1809 break;
1810 case(ContactElementsZero):
1811 out=mesh->ContactElements->Id;
1812 break;
1813 case(ReducedContactElementsZero):
1814 out=mesh->ContactElements->Id;
1815 break;
1816 case(ContactElementsOne):
1817 out=mesh->ContactElements->Id;
1818 break;
1819 case(ReducedContactElementsOne):
1820 out=mesh->ContactElements->Id;
1821 break;
1822 case(DegreesOfFreedom):
1823 out=mesh->Nodes->degreesOfFreedomId;
1824 break;
1825 case(ReducedDegreesOfFreedom):
1826 out=mesh->Nodes->reducedDegreesOfFreedomId;
1827 break;
1828 default:
1829 stringstream temp;
1830 temp << "Error - Invalid function space type: " << functionSpaceType << " for domain: " << getDescription();
1831 throw FinleyAdapterException(temp.str());
1832 break;
1833 }
1834 return out;
1835 }
1836 int MeshAdapter::getTagFromSampleNo(int functionSpaceType, int sampleNo) const
1837 {
1838 int out=0;
1839 Finley_Mesh* mesh=m_finleyMesh.get();
1840 switch (functionSpaceType) {
1841 case(Nodes):
1842 out=mesh->Nodes->Tag[sampleNo];
1843 break;
1844 case(ReducedNodes):
1845 throw FinleyAdapterException(" Error - ReducedNodes does not support tags.");
1846 break;
1847 case(Elements):
1848 out=mesh->Elements->Tag[sampleNo];
1849 break;
1850 case(ReducedElements):
1851 out=mesh->Elements->Tag[sampleNo];
1852 break;
1853 case(FaceElements):
1854 out=mesh->FaceElements->Tag[sampleNo];
1855 break;
1856 case(ReducedFaceElements):
1857 out=mesh->FaceElements->Tag[sampleNo];
1858 break;
1859 case(Points):
1860 out=mesh->Points->Tag[sampleNo];
1861 break;
1862 case(ContactElementsZero):
1863 out=mesh->ContactElements->Tag[sampleNo];
1864 break;
1865 case(ReducedContactElementsZero):
1866 out=mesh->ContactElements->Tag[sampleNo];
1867 break;
1868 case(ContactElementsOne):
1869 out=mesh->ContactElements->Tag[sampleNo];
1870 break;
1871 case(ReducedContactElementsOne):
1872 out=mesh->ContactElements->Tag[sampleNo];
1873 break;
1874 case(DegreesOfFreedom):
1875 throw FinleyAdapterException(" Error - DegreesOfFreedom does not support tags.");
1876 break;
1877 case(ReducedDegreesOfFreedom):
1878 throw FinleyAdapterException(" Error - ReducedDegreesOfFreedom does not support tags.");
1879 break;
1880 default:
1881 stringstream temp;
1882 temp << "Error - Invalid function space type: " << functionSpaceType << " for domain: " << getDescription();
1883 throw FinleyAdapterException(temp.str());
1884 break;
1885 }
1886 return out;
1887 }
1888
1889
1890 void MeshAdapter::setTags(const int functionSpaceType, const int newTag, const escript::Data& mask) const
1891 {
1892 Finley_Mesh* mesh=m_finleyMesh.get();
1893 escriptDataC tmp=mask.getDataC();
1894 switch(functionSpaceType) {
1895 case(Nodes):
1896 Finley_NodeFile_setTags(mesh->Nodes,newTag,&tmp);
1897 break;
1898 case(ReducedNodes):
1899 throw FinleyAdapterException("Error - ReducedNodes does not support tags");
1900 break;
1901 case(DegreesOfFreedom):
1902 throw FinleyAdapterException("Error - DegreesOfFreedom does not support tags");
1903 break;
1904 case(ReducedDegreesOfFreedom):
1905 throw FinleyAdapterException("Error - ReducedDegreesOfFreedom does not support tags");
1906 break;
1907 case(Elements):
1908 Finley_ElementFile_setTags(mesh->Elements,newTag,&tmp);
1909 break;
1910 case(ReducedElements):
1911 Finley_ElementFile_setTags(mesh->Elements,newTag,&tmp);
1912 break;
1913 case(FaceElements):
1914 Finley_ElementFile_setTags(mesh->FaceElements,newTag,&tmp);
1915 break;
1916 case(ReducedFaceElements):
1917 Finley_ElementFile_setTags(mesh->FaceElements,newTag,&tmp);
1918 break;
1919 case(Points):
1920 Finley_ElementFile_setTags(mesh->Points,newTag,&tmp);
1921 break;
1922 case(ContactElementsZero):
1923 Finley_ElementFile_setTags(mesh->ContactElements,newTag,&tmp);
1924 break;
1925 case(ReducedContactElementsZero):
1926 Finley_ElementFile_setTags(mesh->ContactElements,newTag,&tmp);
1927 break;
1928 case(ContactElementsOne):
1929 Finley_ElementFile_setTags(mesh->ContactElements,newTag,&tmp);
1930 break;
1931 case(ReducedContactElementsOne):
1932 Finley_ElementFile_setTags(mesh->ContactElements,newTag,&tmp);
1933 break;
1934 default:
1935 stringstream temp;
1936 temp << "Error - Finley does not know anything about function space type " << functionSpaceType;
1937 throw FinleyAdapterException(temp.str());
1938 }
1939 checkFinleyError();
1940 return;
1941 }
1942
1943 void MeshAdapter::setTagMap(const std::string& name, int tag)
1944 {
1945 Finley_Mesh* mesh=m_finleyMesh.get();
1946 Finley_Mesh_addTagMap(mesh, name.c_str(),tag);
1947 checkPasoError();
1948 // throwStandardException("MeshAdapter::set TagMap is not implemented.");
1949 }
1950
1951 int MeshAdapter::getTag(const std::string& name) const
1952 {
1953 Finley_Mesh* mesh=m_finleyMesh.get();
1954 int tag=0;
1955 tag=Finley_Mesh_getTag(mesh, name.c_str());
1956 checkPasoError();
1957 // throwStandardException("MeshAdapter::getTag is not implemented.");
1958 return tag;
1959 }
1960
1961 bool MeshAdapter::isValidTagName(const std::string& name) const
1962 {
1963 Finley_Mesh* mesh=m_finleyMesh.get();
1964 return Finley_Mesh_isValidTagName(mesh,name.c_str());
1965 }
1966
1967 std::string MeshAdapter::showTagNames() const
1968 {
1969 stringstream temp;
1970 Finley_Mesh* mesh=m_finleyMesh.get();
1971 Finley_TagMap* tag_map=mesh->TagMap;
1972 while (tag_map) {
1973 temp << tag_map->name;
1974 tag_map=tag_map->next;
1975 if (tag_map) temp << ", ";
1976 }
1977 return temp.str();
1978 }
1979
1980 } // end of namespace

Properties

Name Value
svn:eol-style native
svn:keywords Author Date Id Revision

  ViewVC Help
Powered by ViewVC 1.1.26