/[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 1375 - (show annotations)
Wed Jan 9 00:15:05 2008 UTC (11 years, 8 months ago) by gross
File size: 86045 byte(s)
bug in interpolation at reduced face elements fixed.
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 // adds PDE of second order into a transport problem
814 //
815 void MeshAdapter::addPDEToTransportProblem(
816 TransportProblemAdapter& tp, escript::Data& source, const escript::Data& M,
817 const escript::Data& A, const escript::Data& B, const escript::Data& C,const escript::Data& D,const escript::Data& X,const escript::Data& Y,
818 const escript::Data& d, const escript::Data& y,
819 const escript::Data& d_contact,const escript::Data& y_contact) const
820 {
821 DataArrayView::ShapeType shape;
822 source.expand();
823 escript:: Data tmp(0.0,M.getDataPointShape(),tp.getFunctionSpace(),true);
824 escriptDataC _source=source.getDataC();
825 escriptDataC _tmp=tmp.getDataC();
826 escriptDataC _M=M.getDataC();
827 escriptDataC _A=A.getDataC();
828 escriptDataC _B=B.getDataC();
829 escriptDataC _C=C.getDataC();
830 escriptDataC _D=D.getDataC();
831 escriptDataC _X=X.getDataC();
832 escriptDataC _Y=Y.getDataC();
833 escriptDataC _d=d.getDataC();
834 escriptDataC _y=y.getDataC();
835 escriptDataC _d_contact=d_contact.getDataC();
836 escriptDataC _y_contact=y_contact.getDataC();
837
838 Finley_Mesh* mesh=m_finleyMesh.get();
839 Paso_FCTransportProblem* _tp = tp.getPaso_FCTransportProblem();
840
841
842 Finley_Assemble_LumpedSystem(mesh->Nodes,mesh->Elements,&_tmp, &_M);
843 checkFinleyError();
844 /* add mass matix to lumped mass matrix of transport problem */
845 double* tmp_prt=getSampleData(&_tmp,0);
846 int i;
847 int n=Paso_FCTransportProblem_getTotalNumRows(_tp);
848 #pragma omp parallel for private(i) schedule(static)
849 for (i=0;i<n ;++i) _tp->lumped_mass_matrix[i]+=tmp_prt[i];
850
851 Finley_Assemble_PDE(mesh->Nodes,mesh->Elements,_tp->transport_matrix, &_source, &_A, 0, 0, &_D, &_X, &_Y );
852 checkFinleyError();
853
854 Finley_Assemble_PDE(mesh->Nodes,mesh->Elements,_tp->flux_matrix, &_source, 0, &_B, &_C, 0, 0, 0 );
855 checkFinleyError();
856
857 Finley_Assemble_PDE(mesh->Nodes,mesh->FaceElements, _tp->transport_matrix, &_source, 0, 0, 0, &_d, 0, &_y );
858 checkFinleyError();
859
860 Finley_Assemble_PDE(mesh->Nodes,mesh->ContactElements, _tp->transport_matrix, &_source , 0, 0, 0, &_d_contact, 0, &_y_contact );
861 checkFinleyError();
862 }
863
864 //
865 // interpolates data between different function spaces:
866 //
867 void MeshAdapter::interpolateOnDomain(escript::Data& target,const escript::Data& in) const
868 {
869 const MeshAdapter& inDomain=dynamic_cast<const MeshAdapter&>(in.getFunctionSpace().getDomain());
870 const MeshAdapter& targetDomain=dynamic_cast<const MeshAdapter&>(target.getFunctionSpace().getDomain());
871 if (inDomain!=*this)
872 throw FinleyAdapterException("Error - Illegal domain of interpolant.");
873 if (targetDomain!=*this)
874 throw FinleyAdapterException("Error - Illegal domain of interpolation target.");
875
876 Finley_Mesh* mesh=m_finleyMesh.get();
877 escriptDataC _target=target.getDataC();
878 escriptDataC _in=in.getDataC();
879 switch(in.getFunctionSpace().getTypeCode()) {
880 case(Nodes):
881 switch(target.getFunctionSpace().getTypeCode()) {
882 case(Nodes):
883 case(ReducedNodes):
884 case(DegreesOfFreedom):
885 case(ReducedDegreesOfFreedom):
886 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
887 break;
888 case(Elements):
889 case(ReducedElements):
890 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in,&_target);
891 break;
892 case(FaceElements):
893 case(ReducedFaceElements):
894 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in,&_target);
895 break;
896 case(Points):
897 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in,&_target);
898 break;
899 case(ContactElementsZero):
900 case(ReducedContactElementsZero):
901 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
902 break;
903 case(ContactElementsOne):
904 case(ReducedContactElementsOne):
905 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
906 break;
907 default:
908 stringstream temp;
909 temp << "Error - Interpolation on Domain: Finley does not know anything about function space type " << target.getFunctionSpace().getTypeCode();
910 throw FinleyAdapterException(temp.str());
911 break;
912 }
913 break;
914 case(ReducedNodes):
915 switch(target.getFunctionSpace().getTypeCode()) {
916 case(Nodes):
917 case(ReducedNodes):
918 case(DegreesOfFreedom):
919 case(ReducedDegreesOfFreedom):
920 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
921 break;
922 case(Elements):
923 case(ReducedElements):
924 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in,&_target);
925 break;
926 case(FaceElements):
927 case(ReducedFaceElements):
928 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in,&_target);
929 break;
930 case(Points):
931 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in,&_target);
932 break;
933 case(ContactElementsZero):
934 case(ReducedContactElementsZero):
935 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
936 break;
937 case(ContactElementsOne):
938 case(ReducedContactElementsOne):
939 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
940 break;
941 default:
942 stringstream temp;
943 temp << "Error - Interpolation on Domain: Finley does not know anything about function space type " << target.getFunctionSpace().getTypeCode();
944 throw FinleyAdapterException(temp.str());
945 break;
946 }
947 break;
948 case(Elements):
949 if (target.getFunctionSpace().getTypeCode()==Elements) {
950 Finley_Assemble_CopyElementData(mesh->Elements,&_target,&_in);
951 } else if (target.getFunctionSpace().getTypeCode()==ReducedElements) {
952 Finley_Assemble_AverageElementData(mesh->Elements,&_target,&_in);
953 } else {
954 throw FinleyAdapterException("Error - No interpolation with data on elements possible.");
955 }
956 break;
957 case(ReducedElements):
958 if (target.getFunctionSpace().getTypeCode()==ReducedElements) {
959 Finley_Assemble_CopyElementData(mesh->Elements,&_target,&_in);
960 } else {
961 throw FinleyAdapterException("Error - No interpolation with data on elements with reduced integration order possible.");
962 }
963 break;
964 case(FaceElements):
965 if (target.getFunctionSpace().getTypeCode()==FaceElements) {
966 Finley_Assemble_CopyElementData(mesh->FaceElements,&_target,&_in);
967 } else if (target.getFunctionSpace().getTypeCode()==ReducedFaceElements) {
968 Finley_Assemble_AverageElementData(mesh->FaceElements,&_target,&_in);
969 } else {
970 throw FinleyAdapterException("Error - No interpolation with data on face elements possible.");
971 }
972 break;
973 case(ReducedFaceElements):
974 if (target.getFunctionSpace().getTypeCode()==ReducedFaceElements) {
975 Finley_Assemble_CopyElementData(mesh->FaceElements,&_target,&_in);
976 } else {
977 throw FinleyAdapterException("Error - No interpolation with data on face elements with reduced integration order possible.");
978 }
979 break;
980 case(Points):
981 if (target.getFunctionSpace().getTypeCode()==Points) {
982 Finley_Assemble_CopyElementData(mesh->Points,&_target,&_in);
983 } else {
984 throw FinleyAdapterException("Error - No interpolation with data on points possible.");
985 }
986 break;
987 case(ContactElementsZero):
988 case(ContactElementsOne):
989 if (target.getFunctionSpace().getTypeCode()==ContactElementsZero || target.getFunctionSpace().getTypeCode()==ContactElementsOne) {
990 Finley_Assemble_CopyElementData(mesh->ContactElements,&_target,&_in);
991 } else if (target.getFunctionSpace().getTypeCode()==ReducedContactElementsZero || target.getFunctionSpace().getTypeCode()==ReducedContactElementsOne) {
992 Finley_Assemble_AverageElementData(mesh->ContactElements,&_target,&_in);
993 } else {
994 throw FinleyAdapterException("Error - No interpolation with data on contact elements possible.");
995 }
996 break;
997 case(ReducedContactElementsZero):
998 case(ReducedContactElementsOne):
999 if (target.getFunctionSpace().getTypeCode()==ReducedContactElementsZero || target.getFunctionSpace().getTypeCode()==ReducedContactElementsOne) {
1000 Finley_Assemble_CopyElementData(mesh->ContactElements,&_target,&_in);
1001 } else {
1002 throw FinleyAdapterException("Error - No interpolation with data on contact elements with reduced integration order possible.");
1003 }
1004 break;
1005 case(DegreesOfFreedom):
1006 switch(target.getFunctionSpace().getTypeCode()) {
1007 case(ReducedDegreesOfFreedom):
1008 case(DegreesOfFreedom):
1009 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
1010 break;
1011
1012 case(Nodes):
1013 case(ReducedNodes):
1014 if (getMPISize()>1) {
1015 escript::Data temp=escript::Data(in);
1016 temp.expand();
1017 escriptDataC _in2 = temp.getDataC();
1018 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in2);
1019 } else {
1020 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
1021 }
1022 break;
1023 case(Elements):
1024 case(ReducedElements):
1025 if (getMPISize()>1) {
1026 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
1027 escriptDataC _in2 = temp.getDataC();
1028 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in2,&_target);
1029 } else {
1030 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in,&_target);
1031 }
1032 break;
1033 case(FaceElements):
1034 case(ReducedFaceElements):
1035 if (getMPISize()>1) {
1036 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
1037 escriptDataC _in2 = temp.getDataC();
1038 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in2,&_target);
1039
1040 } else {
1041 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in,&_target);
1042 }
1043 break;
1044 case(Points):
1045 if (getMPISize()>1) {
1046 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
1047 escriptDataC _in2 = temp.getDataC();
1048 } else {
1049 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in,&_target);
1050 }
1051 break;
1052 case(ContactElementsZero):
1053 case(ContactElementsOne):
1054 case(ReducedContactElementsZero):
1055 case(ReducedContactElementsOne):
1056 if (getMPISize()>1) {
1057 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
1058 escriptDataC _in2 = temp.getDataC();
1059 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in2,&_target);
1060 } else {
1061 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
1062 }
1063 break;
1064 default:
1065 stringstream temp;
1066 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << target.getFunctionSpace().getTypeCode();
1067 throw FinleyAdapterException(temp.str());
1068 break;
1069 }
1070 break;
1071 case(ReducedDegreesOfFreedom):
1072 switch(target.getFunctionSpace().getTypeCode()) {
1073 case(Nodes):
1074 throw FinleyAdapterException("Error - Finley does not support interpolation from reduced degrees of freedom to mesh nodes.");
1075 break;
1076 case(ReducedNodes):
1077 if (getMPISize()>1) {
1078 escript::Data temp=escript::Data(in);
1079 temp.expand();
1080 escriptDataC _in2 = temp.getDataC();
1081 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in2);
1082 } else {
1083 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
1084 }
1085 break;
1086 case(DegreesOfFreedom):
1087 throw FinleyAdapterException("Error - Finley does not support interpolation from reduced degrees of freedom to degrees of freedom");
1088 break;
1089 case(ReducedDegreesOfFreedom):
1090 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
1091 break;
1092 case(Elements):
1093 case(ReducedElements):
1094 if (getMPISize()>1) {
1095 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
1096 escriptDataC _in2 = temp.getDataC();
1097 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in2,&_target);
1098 } else {
1099 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in,&_target);
1100 }
1101 break;
1102 case(FaceElements):
1103 case(ReducedFaceElements):
1104 if (getMPISize()>1) {
1105 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
1106 escriptDataC _in2 = temp.getDataC();
1107 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in2,&_target);
1108 } else {
1109 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in,&_target);
1110 }
1111 break;
1112 case(Points):
1113 if (getMPISize()>1) {
1114 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
1115 escriptDataC _in2 = temp.getDataC();
1116 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in2,&_target);
1117 } else {
1118 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in,&_target);
1119 }
1120 break;
1121 case(ContactElementsZero):
1122 case(ContactElementsOne):
1123 case(ReducedContactElementsZero):
1124 case(ReducedContactElementsOne):
1125 if (getMPISize()>1) {
1126 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
1127 escriptDataC _in2 = temp.getDataC();
1128 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in2,&_target);
1129 } else {
1130 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
1131 }
1132 break;
1133 default:
1134 stringstream temp;
1135 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << target.getFunctionSpace().getTypeCode();
1136 throw FinleyAdapterException(temp.str());
1137 break;
1138 }
1139 break;
1140 default:
1141 stringstream temp;
1142 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type %d" << in.getFunctionSpace().getTypeCode();
1143 throw FinleyAdapterException(temp.str());
1144 break;
1145 }
1146 checkFinleyError();
1147 }
1148
1149 //
1150 // copies the locations of sample points into x:
1151 //
1152 void MeshAdapter::setToX(escript::Data& arg) const
1153 {
1154 const MeshAdapter& argDomain=dynamic_cast<const MeshAdapter&>(arg.getFunctionSpace().getDomain());
1155 if (argDomain!=*this)
1156 throw FinleyAdapterException("Error - Illegal domain of data point locations");
1157 Finley_Mesh* mesh=m_finleyMesh.get();
1158 // in case of values node coordinates we can do the job directly:
1159 if (arg.getFunctionSpace().getTypeCode()==Nodes) {
1160 escriptDataC _arg=arg.getDataC();
1161 Finley_Assemble_NodeCoordinates(mesh->Nodes,&_arg);
1162 } else {
1163 escript::Data tmp_data=Vector(0.0,continuousFunction(asAbstractContinuousDomain()),true);
1164 escriptDataC _tmp_data=tmp_data.getDataC();
1165 Finley_Assemble_NodeCoordinates(mesh->Nodes,&_tmp_data);
1166 // this is then interpolated onto arg:
1167 interpolateOnDomain(arg,tmp_data);
1168 }
1169 checkFinleyError();
1170 }
1171
1172 //
1173 // return the normal vectors at the location of data points as a Data object:
1174 //
1175 void MeshAdapter::setToNormal(escript::Data& normal) const
1176 {
1177 const MeshAdapter& normalDomain=dynamic_cast<const MeshAdapter&>(normal.getFunctionSpace().getDomain());
1178 if (normalDomain!=*this)
1179 throw FinleyAdapterException("Error - Illegal domain of normal locations");
1180 Finley_Mesh* mesh=m_finleyMesh.get();
1181 escriptDataC _normal=normal.getDataC();
1182 switch(normal.getFunctionSpace().getTypeCode()) {
1183 case(Nodes):
1184 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for nodes");
1185 break;
1186 case(ReducedNodes):
1187 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for reduced nodes");
1188 break;
1189 case(Elements):
1190 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for elements");
1191 break;
1192 case(ReducedElements):
1193 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for elements with reduced integration order");
1194 break;
1195 case (FaceElements):
1196 Finley_Assemble_setNormal(mesh->Nodes,mesh->FaceElements,&_normal);
1197 break;
1198 case (ReducedFaceElements):
1199 Finley_Assemble_setNormal(mesh->Nodes,mesh->FaceElements,&_normal);
1200 break;
1201 case(Points):
1202 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for point elements");
1203 break;
1204 case (ContactElementsOne):
1205 case (ContactElementsZero):
1206 Finley_Assemble_setNormal(mesh->Nodes,mesh->ContactElements,&_normal);
1207 break;
1208 case (ReducedContactElementsOne):
1209 case (ReducedContactElementsZero):
1210 Finley_Assemble_setNormal(mesh->Nodes,mesh->ContactElements,&_normal);
1211 break;
1212 case(DegreesOfFreedom):
1213 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for degrees of freedom.");
1214 break;
1215 case(ReducedDegreesOfFreedom):
1216 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for reduced degrees of freedom.");
1217 break;
1218 default:
1219 stringstream temp;
1220 temp << "Error - Normal Vectors: Finley does not know anything about function space type " << normal.getFunctionSpace().getTypeCode();
1221 throw FinleyAdapterException(temp.str());
1222 break;
1223 }
1224 checkFinleyError();
1225 }
1226
1227 //
1228 // interpolates data to other domain:
1229 //
1230 void MeshAdapter::interpolateACross(escript::Data& target,const escript::Data& source) const
1231 {
1232 const MeshAdapter& targetDomain=dynamic_cast<const MeshAdapter&>(target.getFunctionSpace().getDomain());
1233 if (targetDomain!=*this)
1234 throw FinleyAdapterException("Error - Illegal domain of interpolation target");
1235
1236 throw FinleyAdapterException("Error - Finley does not allow interpolation across domains yet.");
1237 }
1238
1239 //
1240 // calculates the integral of a function defined of arg:
1241 //
1242 void MeshAdapter::setToIntegrals(std::vector<double>& integrals,const escript::Data& arg) const
1243 {
1244 const MeshAdapter& argDomain=dynamic_cast<const MeshAdapter&>(arg.getFunctionSpace().getDomain());
1245 if (argDomain!=*this)
1246 throw FinleyAdapterException("Error - Illegal domain of integration kernel");
1247
1248 double blocktimer_start = blocktimer_time();
1249 Finley_Mesh* mesh=m_finleyMesh.get();
1250 escriptDataC _temp;
1251 escript::Data temp;
1252 escriptDataC _arg=arg.getDataC();
1253 switch(arg.getFunctionSpace().getTypeCode()) {
1254 case(Nodes):
1255 temp=escript::Data( arg, function(asAbstractContinuousDomain()) );
1256 _temp=temp.getDataC();
1257 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_temp,&integrals[0]);
1258 break;
1259 case(ReducedNodes):
1260 temp=escript::Data( arg, function(asAbstractContinuousDomain()) );
1261 _temp=temp.getDataC();
1262 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_temp,&integrals[0]);
1263 break;
1264 case(Elements):
1265 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_arg,&integrals[0]);
1266 break;
1267 case(ReducedElements):
1268 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_arg,&integrals[0]);
1269 break;
1270 case(FaceElements):
1271 Finley_Assemble_integrate(mesh->Nodes,mesh->FaceElements,&_arg,&integrals[0]);
1272 break;
1273 case(ReducedFaceElements):
1274 Finley_Assemble_integrate(mesh->Nodes,mesh->FaceElements,&_arg,&integrals[0]);
1275 break;
1276 case(Points):
1277 throw FinleyAdapterException("Error - Integral of data on points is not supported.");
1278 break;
1279 case(ContactElementsZero):
1280 Finley_Assemble_integrate(mesh->Nodes,mesh->ContactElements,&_arg,&integrals[0]);
1281 break;
1282 case(ReducedContactElementsZero):
1283 Finley_Assemble_integrate(mesh->Nodes,mesh->ContactElements,&_arg,&integrals[0]);
1284 break;
1285 case(ContactElementsOne):
1286 Finley_Assemble_integrate(mesh->Nodes,mesh->ContactElements,&_arg,&integrals[0]);
1287 break;
1288 case(ReducedContactElementsOne):
1289 Finley_Assemble_integrate(mesh->Nodes,mesh->ContactElements,&_arg,&integrals[0]);
1290 break;
1291 case(DegreesOfFreedom):
1292 temp=escript::Data( arg, function(asAbstractContinuousDomain()) );
1293 _temp=temp.getDataC();
1294 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_temp,&integrals[0]);
1295 break;
1296 case(ReducedDegreesOfFreedom):
1297 temp=escript::Data( arg, function(asAbstractContinuousDomain()) );
1298 _temp=temp.getDataC();
1299 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_temp,&integrals[0]);
1300 break;
1301 default:
1302 stringstream temp;
1303 temp << "Error - Integrals: Finley does not know anything about function space type " << arg.getFunctionSpace().getTypeCode();
1304 throw FinleyAdapterException(temp.str());
1305 break;
1306 }
1307 checkFinleyError();
1308 blocktimer_increment("integrate()", blocktimer_start);
1309 }
1310
1311 //
1312 // calculates the gradient of arg:
1313 //
1314 void MeshAdapter::setToGradient(escript::Data& grad,const escript::Data& arg) const
1315 {
1316 const MeshAdapter& argDomain=dynamic_cast<const MeshAdapter&>(arg.getFunctionSpace().getDomain());
1317 if (argDomain!=*this)
1318 throw FinleyAdapterException("Error - Illegal domain of gradient argument");
1319 const MeshAdapter& gradDomain=dynamic_cast<const MeshAdapter&>(grad.getFunctionSpace().getDomain());
1320 if (gradDomain!=*this)
1321 throw FinleyAdapterException("Error - Illegal domain of gradient");
1322
1323 Finley_Mesh* mesh=m_finleyMesh.get();
1324 escriptDataC _grad=grad.getDataC();
1325 escriptDataC nodeDataC;
1326 escript::Data temp;
1327 if (getMPISize()>1) {
1328 if( arg.getFunctionSpace().getTypeCode() == DegreesOfFreedom ) {
1329 temp=escript::Data( arg, continuousFunction(asAbstractContinuousDomain()) );
1330 nodeDataC = temp.getDataC();
1331 } else if( arg.getFunctionSpace().getTypeCode() == ReducedDegreesOfFreedom ) {
1332 temp=escript::Data( arg, reducedContinuousFunction(asAbstractContinuousDomain()) );
1333 nodeDataC = temp.getDataC();
1334 } else {
1335 nodeDataC = arg.getDataC();
1336 }
1337 } else {
1338 nodeDataC = arg.getDataC();
1339 }
1340 switch(grad.getFunctionSpace().getTypeCode()) {
1341 case(Nodes):
1342 throw FinleyAdapterException("Error - Gradient at nodes is not supported.");
1343 break;
1344 case(ReducedNodes):
1345 throw FinleyAdapterException("Error - Gradient at reduced nodes is not supported.");
1346 break;
1347 case(Elements):
1348 Finley_Assemble_gradient(mesh->Nodes,mesh->Elements,&_grad,&nodeDataC);
1349 break;
1350 case(ReducedElements):
1351 Finley_Assemble_gradient(mesh->Nodes,mesh->Elements,&_grad,&nodeDataC);
1352 break;
1353 case(FaceElements):
1354 Finley_Assemble_gradient(mesh->Nodes,mesh->FaceElements,&_grad,&nodeDataC);
1355 break;
1356 case(ReducedFaceElements):
1357 Finley_Assemble_gradient(mesh->Nodes,mesh->FaceElements,&_grad,&nodeDataC);
1358 break;
1359 case(Points):
1360 throw FinleyAdapterException("Error - Gradient at points is not supported.");
1361 break;
1362 case(ContactElementsZero):
1363 Finley_Assemble_gradient(mesh->Nodes,mesh->ContactElements,&_grad,&nodeDataC);
1364 break;
1365 case(ReducedContactElementsZero):
1366 Finley_Assemble_gradient(mesh->Nodes,mesh->ContactElements,&_grad,&nodeDataC);
1367 break;
1368 case(ContactElementsOne):
1369 Finley_Assemble_gradient(mesh->Nodes,mesh->ContactElements,&_grad,&nodeDataC);
1370 break;
1371 case(ReducedContactElementsOne):
1372 Finley_Assemble_gradient(mesh->Nodes,mesh->ContactElements,&_grad,&nodeDataC);
1373 break;
1374 case(DegreesOfFreedom):
1375 throw FinleyAdapterException("Error - Gradient at degrees of freedom is not supported.");
1376 break;
1377 case(ReducedDegreesOfFreedom):
1378 throw FinleyAdapterException("Error - Gradient at reduced degrees of freedom is not supported.");
1379 break;
1380 default:
1381 stringstream temp;
1382 temp << "Error - Gradient: Finley does not know anything about function space type " << arg.getFunctionSpace().getTypeCode();
1383 throw FinleyAdapterException(temp.str());
1384 break;
1385 }
1386 checkFinleyError();
1387 }
1388
1389 //
1390 // returns the size of elements:
1391 //
1392 void MeshAdapter::setToSize(escript::Data& size) const
1393 {
1394 Finley_Mesh* mesh=m_finleyMesh.get();
1395 escriptDataC tmp=size.getDataC();
1396 switch(size.getFunctionSpace().getTypeCode()) {
1397 case(Nodes):
1398 throw FinleyAdapterException("Error - Size of nodes is not supported.");
1399 break;
1400 case(ReducedNodes):
1401 throw FinleyAdapterException("Error - Size of reduced nodes is not supported.");
1402 break;
1403 case(Elements):
1404 Finley_Assemble_getSize(mesh->Nodes,mesh->Elements,&tmp);
1405 break;
1406 case(ReducedElements):
1407 Finley_Assemble_getSize(mesh->Nodes,mesh->Elements,&tmp);
1408 break;
1409 case(FaceElements):
1410 Finley_Assemble_getSize(mesh->Nodes,mesh->FaceElements,&tmp);
1411 break;
1412 case(ReducedFaceElements):
1413 Finley_Assemble_getSize(mesh->Nodes,mesh->FaceElements,&tmp);
1414 break;
1415 case(Points):
1416 throw FinleyAdapterException("Error - Size of point elements is not supported.");
1417 break;
1418 case(ContactElementsZero):
1419 case(ContactElementsOne):
1420 Finley_Assemble_getSize(mesh->Nodes,mesh->ContactElements,&tmp);
1421 break;
1422 case(ReducedContactElementsZero):
1423 case(ReducedContactElementsOne):
1424 Finley_Assemble_getSize(mesh->Nodes,mesh->ContactElements,&tmp);
1425 break;
1426 case(DegreesOfFreedom):
1427 throw FinleyAdapterException("Error - Size of degrees of freedom is not supported.");
1428 break;
1429 case(ReducedDegreesOfFreedom):
1430 throw FinleyAdapterException("Error - Size of reduced degrees of freedom is not supported.");
1431 break;
1432 default:
1433 stringstream temp;
1434 temp << "Error - Element size: Finley does not know anything about function space type " << size.getFunctionSpace().getTypeCode();
1435 throw FinleyAdapterException(temp.str());
1436 break;
1437 }
1438 checkFinleyError();
1439 }
1440
1441 // sets the location of nodes:
1442 void MeshAdapter::setNewX(const escript::Data& new_x)
1443 {
1444 Finley_Mesh* mesh=m_finleyMesh.get();
1445 escriptDataC tmp;
1446 const MeshAdapter& newDomain=dynamic_cast<const MeshAdapter&>(new_x.getFunctionSpace().getDomain());
1447 if (newDomain!=*this)
1448 throw FinleyAdapterException("Error - Illegal domain of new point locations");
1449 tmp = new_x.getDataC();
1450 Finley_Mesh_setCoordinates(mesh,&tmp);
1451 checkFinleyError();
1452 }
1453
1454 // saves a data array in openDX format:
1455 void MeshAdapter::saveDX(const std::string& filename,const boost::python::dict& arg) const
1456 {
1457 int MAX_namelength=256;
1458 const int num_data=boost::python::extract<int>(arg.attr("__len__")());
1459 /* win32 refactor */
1460 char* *names = (num_data>0) ? TMPMEMALLOC(num_data,char*) : (char**)NULL;
1461 for(int i=0;i<num_data;i++)
1462 {
1463 names[i] = (MAX_namelength>0) ? TMPMEMALLOC(MAX_namelength,char) : (char*)NULL;
1464 }
1465
1466 char* *c_names = (num_data>0) ? TMPMEMALLOC(num_data,char*) : (char**)NULL;
1467 escriptDataC *data = (num_data>0) ? TMPMEMALLOC(num_data,escriptDataC) : (escriptDataC*)NULL;
1468 escriptDataC* *ptr_data = (num_data>0) ? TMPMEMALLOC(num_data,escriptDataC*) : (escriptDataC**)NULL;
1469
1470 boost::python::list keys=arg.keys();
1471 for (int i=0;i<num_data;++i) {
1472 std::string n=boost::python::extract<std::string>(keys[i]);
1473 escript::Data& d=boost::python::extract<escript::Data&>(arg[keys[i]]);
1474 if (dynamic_cast<const MeshAdapter&>(d.getFunctionSpace().getDomain()) !=*this)
1475 throw FinleyAdapterException("Error in saveVTK: Data must be defined on same Domain");
1476 data[i]=d.getDataC();
1477 ptr_data[i]=&(data[i]);
1478 c_names[i]=&(names[i][0]);
1479 if (n.length()>MAX_namelength-1) {
1480 strncpy(c_names[i],n.c_str(),MAX_namelength-1);
1481 c_names[i][MAX_namelength-1]='\0';
1482 } else {
1483 strcpy(c_names[i],n.c_str());
1484 }
1485 }
1486 Finley_Mesh_saveDX(filename.c_str(),m_finleyMesh.get(),num_data,c_names,ptr_data);
1487 checkFinleyError();
1488
1489 /* win32 refactor */
1490 TMPMEMFREE(c_names);
1491 TMPMEMFREE(data);
1492 TMPMEMFREE(ptr_data);
1493 for(int i=0;i<num_data;i++)
1494 {
1495 TMPMEMFREE(names[i]);
1496 }
1497 TMPMEMFREE(names);
1498
1499 return;
1500 }
1501
1502 // saves a data array in openVTK format:
1503 void MeshAdapter::saveVTK(const std::string& filename,const boost::python::dict& arg) const
1504 {
1505 int MAX_namelength=256;
1506 const int num_data=boost::python::extract<int>(arg.attr("__len__")());
1507 /* win32 refactor */
1508 char* *names = (num_data>0) ? TMPMEMALLOC(num_data,char*) : (char**)NULL;
1509 for(int i=0;i<num_data;i++)
1510 {
1511 names[i] = (MAX_namelength>0) ? TMPMEMALLOC(MAX_namelength,char) : (char*)NULL;
1512 }
1513
1514 char* *c_names = (num_data>0) ? TMPMEMALLOC(num_data,char*) : (char**)NULL;
1515 escriptDataC *data = (num_data>0) ? TMPMEMALLOC(num_data,escriptDataC) : (escriptDataC*)NULL;
1516 escriptDataC* *ptr_data = (num_data>0) ? TMPMEMALLOC(num_data,escriptDataC*) : (escriptDataC**)NULL;
1517
1518 boost::python::list keys=arg.keys();
1519 for (int i=0;i<num_data;++i) {
1520 std::string n=boost::python::extract<std::string>(keys[i]);
1521 escript::Data& d=boost::python::extract<escript::Data&>(arg[keys[i]]);
1522 if (dynamic_cast<const MeshAdapter&>(d.getFunctionSpace().getDomain()) !=*this)
1523 throw FinleyAdapterException("Error in saveVTK: Data must be defined on same Domain");
1524 data[i]=d.getDataC();
1525 ptr_data[i]=&(data[i]);
1526 c_names[i]=&(names[i][0]);
1527 if (n.length()>MAX_namelength-1) {
1528 strncpy(c_names[i],n.c_str(),MAX_namelength-1);
1529 c_names[i][MAX_namelength-1]='\0';
1530 } else {
1531 strcpy(c_names[i],n.c_str());
1532 }
1533 }
1534 Finley_Mesh_saveVTK(filename.c_str(),m_finleyMesh.get(),num_data,c_names,ptr_data);
1535
1536 checkFinleyError();
1537 /* win32 refactor */
1538 TMPMEMFREE(c_names);
1539 TMPMEMFREE(data);
1540 TMPMEMFREE(ptr_data);
1541 for(int i=0;i<num_data;i++)
1542 {
1543 TMPMEMFREE(names[i]);
1544 }
1545 TMPMEMFREE(names);
1546
1547 return;
1548 }
1549
1550
1551 // creates a SystemMatrixAdapter stiffness matrix an initializes it with zeros:
1552 SystemMatrixAdapter MeshAdapter::newSystemMatrix(
1553 const int row_blocksize,
1554 const escript::FunctionSpace& row_functionspace,
1555 const int column_blocksize,
1556 const escript::FunctionSpace& column_functionspace,
1557 const int type) const
1558 {
1559 int reduceRowOrder=0;
1560 int reduceColOrder=0;
1561 // is the domain right?
1562 const MeshAdapter& row_domain=dynamic_cast<const MeshAdapter&>(row_functionspace.getDomain());
1563 if (row_domain!=*this)
1564 throw FinleyAdapterException("Error - domain of row function space does not match the domain of matrix generator.");
1565 const MeshAdapter& col_domain=dynamic_cast<const MeshAdapter&>(column_functionspace.getDomain());
1566 if (col_domain!=*this)
1567 throw FinleyAdapterException("Error - domain of columnn function space does not match the domain of matrix generator.");
1568 // is the function space type right
1569 if (row_functionspace.getTypeCode()==DegreesOfFreedom) {
1570 reduceRowOrder=0;
1571 } else if (row_functionspace.getTypeCode()==ReducedDegreesOfFreedom) {
1572 reduceRowOrder=1;
1573 } else {
1574 throw FinleyAdapterException("Error - illegal function space type for system matrix rows.");
1575 }
1576 if (column_functionspace.getTypeCode()==DegreesOfFreedom) {
1577 reduceColOrder=0;
1578 } else if (column_functionspace.getTypeCode()==ReducedDegreesOfFreedom) {
1579 reduceColOrder=1;
1580 } else {
1581 throw FinleyAdapterException("Error - illegal function space type for system matrix columns.");
1582 }
1583 // generate matrix:
1584
1585 Paso_SystemMatrixPattern* fsystemMatrixPattern=Finley_getPattern(getFinley_Mesh(),reduceRowOrder,reduceColOrder);
1586 checkFinleyError();
1587 Paso_SystemMatrix* fsystemMatrix;
1588 int trilinos = 0;
1589 if (trilinos) {
1590 #ifdef TRILINOS
1591 /* Allocation Epetra_VrbMatrix here */
1592 #endif
1593 }
1594 else {
1595 fsystemMatrix=Paso_SystemMatrix_alloc(type,fsystemMatrixPattern,row_blocksize,column_blocksize);
1596 }
1597 checkPasoError();
1598 Paso_SystemMatrixPattern_free(fsystemMatrixPattern);
1599 return SystemMatrixAdapter(fsystemMatrix,row_blocksize,row_functionspace,column_blocksize,column_functionspace);
1600 }
1601 // creates a TransportProblemAdapter
1602 TransportProblemAdapter MeshAdapter::newTransportProblem(
1603 const double theta,
1604 const double dt_max,
1605 const int blocksize,
1606 const escript::FunctionSpace& functionspace,
1607 const int type) const
1608 {
1609 int reduceOrder=0;
1610 // is the domain right?
1611 const MeshAdapter& domain=dynamic_cast<const MeshAdapter&>(functionspace.getDomain());
1612 if (domain!=*this)
1613 throw FinleyAdapterException("Error - domain of function space does not match the domain of transport problem generator.");
1614 // is the function space type right
1615 if (functionspace.getTypeCode()==DegreesOfFreedom) {
1616 reduceOrder=0;
1617 } else if (functionspace.getTypeCode()==ReducedDegreesOfFreedom) {
1618 reduceOrder=1;
1619 } else {
1620 throw FinleyAdapterException("Error - illegal function space type for system matrix rows.");
1621 }
1622 // generate matrix:
1623
1624 Paso_SystemMatrixPattern* fsystemMatrixPattern=Finley_getPattern(getFinley_Mesh(),reduceOrder,reduceOrder);
1625 checkFinleyError();
1626 Paso_FCTransportProblem* transportProblem;
1627 transportProblem=Paso_FCTransportProblem_alloc(theta,dt_max,fsystemMatrixPattern,blocksize);
1628 checkPasoError();
1629 Paso_SystemMatrixPattern_free(fsystemMatrixPattern);
1630 return TransportProblemAdapter(transportProblem,theta,dt_max,blocksize,functionspace);
1631 }
1632
1633 //
1634 // vtkObject MeshAdapter::createVtkObject() const
1635 // TODO:
1636 //
1637
1638 //
1639 // returns true if data at the atom_type is considered as being cell centered:
1640 bool MeshAdapter::isCellOriented(int functionSpaceCode) const
1641 {
1642 switch(functionSpaceCode) {
1643 case(Nodes):
1644 case(DegreesOfFreedom):
1645 case(ReducedDegreesOfFreedom):
1646 return false;
1647 break;
1648 case(Elements):
1649 case(FaceElements):
1650 case(Points):
1651 case(ContactElementsZero):
1652 case(ContactElementsOne):
1653 case(ReducedElements):
1654 case(ReducedFaceElements):
1655 case(ReducedContactElementsZero):
1656 case(ReducedContactElementsOne):
1657 return true;
1658 break;
1659 default:
1660 stringstream temp;
1661 temp << "Error - Cell: Finley does not know anything about function space type " << functionSpaceCode;
1662 throw FinleyAdapterException(temp.str());
1663 break;
1664 }
1665 checkFinleyError();
1666 return false;
1667 }
1668
1669 bool MeshAdapter::probeInterpolationOnDomain(int functionSpaceType_source,int functionSpaceType_target) const
1670 {
1671 switch(functionSpaceType_source) {
1672 case(Nodes):
1673 switch(functionSpaceType_target) {
1674 case(Nodes):
1675 case(ReducedNodes):
1676 case(ReducedDegreesOfFreedom):
1677 case(DegreesOfFreedom):
1678 case(Elements):
1679 case(ReducedElements):
1680 case(FaceElements):
1681 case(ReducedFaceElements):
1682 case(Points):
1683 case(ContactElementsZero):
1684 case(ReducedContactElementsZero):
1685 case(ContactElementsOne):
1686 case(ReducedContactElementsOne):
1687 return true;
1688 default:
1689 stringstream temp;
1690 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_target;
1691 throw FinleyAdapterException(temp.str());
1692 }
1693 break;
1694 case(ReducedNodes):
1695 switch(functionSpaceType_target) {
1696 case(ReducedNodes):
1697 case(ReducedDegreesOfFreedom):
1698 case(Elements):
1699 case(ReducedElements):
1700 case(FaceElements):
1701 case(ReducedFaceElements):
1702 case(Points):
1703 case(ContactElementsZero):
1704 case(ReducedContactElementsZero):
1705 case(ContactElementsOne):
1706 case(ReducedContactElementsOne):
1707 return true;
1708 case(Nodes):
1709 case(DegreesOfFreedom):
1710 return false;
1711 default:
1712 stringstream temp;
1713 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_target;
1714 throw FinleyAdapterException(temp.str());
1715 }
1716 break;
1717 case(Elements):
1718 if (functionSpaceType_target==Elements) {
1719 return true;
1720 } else if (functionSpaceType_target==ReducedElements) {
1721 return true;
1722 } else {
1723 return false;
1724 }
1725 case(ReducedElements):
1726 if (functionSpaceType_target==ReducedElements) {
1727 return true;
1728 } else {
1729 return false;
1730 }
1731 case(FaceElements):
1732 if (functionSpaceType_target==FaceElements) {
1733 return true;
1734 } else if (functionSpaceType_target==ReducedFaceElements) {
1735 return true;
1736 } else {
1737 return false;
1738 }
1739 case(ReducedFaceElements):
1740 if (functionSpaceType_target==ReducedFaceElements) {
1741 return true;
1742 } else {
1743 return false;
1744 }
1745 case(Points):
1746 if (functionSpaceType_target==Points) {
1747 return true;
1748 } else {
1749 return false;
1750 }
1751 case(ContactElementsZero):
1752 case(ContactElementsOne):
1753 if (functionSpaceType_target==ContactElementsZero || functionSpaceType_target==ContactElementsOne) {
1754 return true;
1755 } else if (functionSpaceType_target==ReducedContactElementsZero || functionSpaceType_target==ReducedContactElementsOne) {
1756 return true;
1757 } else {
1758 return false;
1759 }
1760 case(ReducedContactElementsZero):
1761 case(ReducedContactElementsOne):
1762 if (functionSpaceType_target==ReducedContactElementsZero || functionSpaceType_target==ReducedContactElementsOne) {
1763 return true;
1764 } else {
1765 return false;
1766 }
1767 case(DegreesOfFreedom):
1768 switch(functionSpaceType_target) {
1769 case(ReducedDegreesOfFreedom):
1770 case(DegreesOfFreedom):
1771 case(Nodes):
1772 case(ReducedNodes):
1773 case(Elements):
1774 case(ReducedElements):
1775 case(Points):
1776 case(FaceElements):
1777 case(ReducedFaceElements):
1778 case(ContactElementsZero):
1779 case(ReducedContactElementsZero):
1780 case(ContactElementsOne):
1781 case(ReducedContactElementsOne):
1782 return true;
1783 default:
1784 stringstream temp;
1785 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_target;
1786 throw FinleyAdapterException(temp.str());
1787 }
1788 break;
1789 case(ReducedDegreesOfFreedom):
1790 switch(functionSpaceType_target) {
1791 case(ReducedDegreesOfFreedom):
1792 case(ReducedNodes):
1793 case(Elements):
1794 case(ReducedElements):
1795 case(FaceElements):
1796 case(ReducedFaceElements):
1797 case(Points):
1798 case(ContactElementsZero):
1799 case(ReducedContactElementsZero):
1800 case(ContactElementsOne):
1801 case(ReducedContactElementsOne):
1802 return true;
1803 case(Nodes):
1804 case(DegreesOfFreedom):
1805 return false;
1806 default:
1807 stringstream temp;
1808 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_target;
1809 throw FinleyAdapterException(temp.str());
1810 }
1811 break;
1812 default:
1813 stringstream temp;
1814 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_source;
1815 throw FinleyAdapterException(temp.str());
1816 break;
1817 }
1818 checkFinleyError();
1819 return false;
1820 }
1821
1822 bool MeshAdapter::probeInterpolationACross(int functionSpaceType_source,const AbstractDomain& targetDomain, int functionSpaceType_target) const
1823 {
1824 return false;
1825 }
1826
1827 bool MeshAdapter::operator==(const AbstractDomain& other) const
1828 {
1829 const MeshAdapter* temp=dynamic_cast<const MeshAdapter*>(&other);
1830 if (temp!=0) {
1831 return (m_finleyMesh==temp->m_finleyMesh);
1832 } else {
1833 return false;
1834 }
1835 }
1836
1837 bool MeshAdapter::operator!=(const AbstractDomain& other) const
1838 {
1839 return !(operator==(other));
1840 }
1841
1842 int MeshAdapter::getSystemMatrixTypeId(const int solver, const int package, const bool symmetry) const
1843 {
1844 int out=Paso_SystemMatrix_getSystemMatrixTypeId(SystemMatrixAdapter::mapOptionToPaso(solver),SystemMatrixAdapter::mapOptionToPaso(package),symmetry?1:0);
1845 checkPasoError();
1846 return out;
1847 }
1848
1849 escript::Data MeshAdapter::getX() const
1850 {
1851 return continuousFunction(asAbstractContinuousDomain()).getX();
1852 }
1853
1854 escript::Data MeshAdapter::getNormal() const
1855 {
1856 return functionOnBoundary(asAbstractContinuousDomain()).getNormal();
1857 }
1858
1859 escript::Data MeshAdapter::getSize() const
1860 {
1861 return function(asAbstractContinuousDomain()).getSize();
1862 }
1863
1864 int* MeshAdapter::borrowSampleReferenceIDs(int functionSpaceType) const
1865 {
1866 int *out=0,i;
1867 Finley_Mesh* mesh=m_finleyMesh.get();
1868 switch (functionSpaceType) {
1869 case(Nodes):
1870 out=mesh->Nodes->Id;
1871 break;
1872 case(ReducedNodes):
1873 out=mesh->Nodes->reducedNodesId;
1874 break;
1875 case(Elements):
1876 out=mesh->Elements->Id;
1877 break;
1878 case(ReducedElements):
1879 out=mesh->Elements->Id;
1880 break;
1881 case(FaceElements):
1882 out=mesh->FaceElements->Id;
1883 break;
1884 case(ReducedFaceElements):
1885 out=mesh->FaceElements->Id;
1886 break;
1887 case(Points):
1888 out=mesh->Points->Id;
1889 break;
1890 case(ContactElementsZero):
1891 out=mesh->ContactElements->Id;
1892 break;
1893 case(ReducedContactElementsZero):
1894 out=mesh->ContactElements->Id;
1895 break;
1896 case(ContactElementsOne):
1897 out=mesh->ContactElements->Id;
1898 break;
1899 case(ReducedContactElementsOne):
1900 out=mesh->ContactElements->Id;
1901 break;
1902 case(DegreesOfFreedom):
1903 out=mesh->Nodes->degreesOfFreedomId;
1904 break;
1905 case(ReducedDegreesOfFreedom):
1906 out=mesh->Nodes->reducedDegreesOfFreedomId;
1907 break;
1908 default:
1909 stringstream temp;
1910 temp << "Error - Invalid function space type: " << functionSpaceType << " for domain: " << getDescription();
1911 throw FinleyAdapterException(temp.str());
1912 break;
1913 }
1914 return out;
1915 }
1916 int MeshAdapter::getTagFromSampleNo(int functionSpaceType, int sampleNo) const
1917 {
1918 int out=0;
1919 Finley_Mesh* mesh=m_finleyMesh.get();
1920 switch (functionSpaceType) {
1921 case(Nodes):
1922 out=mesh->Nodes->Tag[sampleNo];
1923 break;
1924 case(ReducedNodes):
1925 throw FinleyAdapterException(" Error - ReducedNodes does not support tags.");
1926 break;
1927 case(Elements):
1928 out=mesh->Elements->Tag[sampleNo];
1929 break;
1930 case(ReducedElements):
1931 out=mesh->Elements->Tag[sampleNo];
1932 break;
1933 case(FaceElements):
1934 out=mesh->FaceElements->Tag[sampleNo];
1935 break;
1936 case(ReducedFaceElements):
1937 out=mesh->FaceElements->Tag[sampleNo];
1938 break;
1939 case(Points):
1940 out=mesh->Points->Tag[sampleNo];
1941 break;
1942 case(ContactElementsZero):
1943 out=mesh->ContactElements->Tag[sampleNo];
1944 break;
1945 case(ReducedContactElementsZero):
1946 out=mesh->ContactElements->Tag[sampleNo];
1947 break;
1948 case(ContactElementsOne):
1949 out=mesh->ContactElements->Tag[sampleNo];
1950 break;
1951 case(ReducedContactElementsOne):
1952 out=mesh->ContactElements->Tag[sampleNo];
1953 break;
1954 case(DegreesOfFreedom):
1955 throw FinleyAdapterException(" Error - DegreesOfFreedom does not support tags.");
1956 break;
1957 case(ReducedDegreesOfFreedom):
1958 throw FinleyAdapterException(" Error - ReducedDegreesOfFreedom does not support tags.");
1959 break;
1960 default:
1961 stringstream temp;
1962 temp << "Error - Invalid function space type: " << functionSpaceType << " for domain: " << getDescription();
1963 throw FinleyAdapterException(temp.str());
1964 break;
1965 }
1966 return out;
1967 }
1968
1969
1970 void MeshAdapter::setTags(const int functionSpaceType, const int newTag, const escript::Data& mask) const
1971 {
1972 Finley_Mesh* mesh=m_finleyMesh.get();
1973 escriptDataC tmp=mask.getDataC();
1974 switch(functionSpaceType) {
1975 case(Nodes):
1976 Finley_NodeFile_setTags(mesh->Nodes,newTag,&tmp);
1977 break;
1978 case(ReducedNodes):
1979 throw FinleyAdapterException("Error - ReducedNodes does not support tags");
1980 break;
1981 case(DegreesOfFreedom):
1982 throw FinleyAdapterException("Error - DegreesOfFreedom does not support tags");
1983 break;
1984 case(ReducedDegreesOfFreedom):
1985 throw FinleyAdapterException("Error - ReducedDegreesOfFreedom does not support tags");
1986 break;
1987 case(Elements):
1988 Finley_ElementFile_setTags(mesh->Elements,newTag,&tmp);
1989 break;
1990 case(ReducedElements):
1991 Finley_ElementFile_setTags(mesh->Elements,newTag,&tmp);
1992 break;
1993 case(FaceElements):
1994 Finley_ElementFile_setTags(mesh->FaceElements,newTag,&tmp);
1995 break;
1996 case(ReducedFaceElements):
1997 Finley_ElementFile_setTags(mesh->FaceElements,newTag,&tmp);
1998 break;
1999 case(Points):
2000 Finley_ElementFile_setTags(mesh->Points,newTag,&tmp);
2001 break;
2002 case(ContactElementsZero):
2003 Finley_ElementFile_setTags(mesh->ContactElements,newTag,&tmp);
2004 break;
2005 case(ReducedContactElementsZero):
2006 Finley_ElementFile_setTags(mesh->ContactElements,newTag,&tmp);
2007 break;
2008 case(ContactElementsOne):
2009 Finley_ElementFile_setTags(mesh->ContactElements,newTag,&tmp);
2010 break;
2011 case(ReducedContactElementsOne):
2012 Finley_ElementFile_setTags(mesh->ContactElements,newTag,&tmp);
2013 break;
2014 default:
2015 stringstream temp;
2016 temp << "Error - Finley does not know anything about function space type " << functionSpaceType;
2017 throw FinleyAdapterException(temp.str());
2018 }
2019 checkFinleyError();
2020 return;
2021 }
2022
2023 void MeshAdapter::setTagMap(const std::string& name, int tag)
2024 {
2025 Finley_Mesh* mesh=m_finleyMesh.get();
2026 Finley_Mesh_addTagMap(mesh, name.c_str(),tag);
2027 checkPasoError();
2028 // throwStandardException("MeshAdapter::set TagMap is not implemented.");
2029 }
2030
2031 int MeshAdapter::getTag(const std::string& name) const
2032 {
2033 Finley_Mesh* mesh=m_finleyMesh.get();
2034 int tag=0;
2035 tag=Finley_Mesh_getTag(mesh, name.c_str());
2036 checkPasoError();
2037 // throwStandardException("MeshAdapter::getTag is not implemented.");
2038 return tag;
2039 }
2040
2041 bool MeshAdapter::isValidTagName(const std::string& name) const
2042 {
2043 Finley_Mesh* mesh=m_finleyMesh.get();
2044 return Finley_Mesh_isValidTagName(mesh,name.c_str());
2045 }
2046
2047 std::string MeshAdapter::showTagNames() const
2048 {
2049 stringstream temp;
2050 Finley_Mesh* mesh=m_finleyMesh.get();
2051 Finley_TagMap* tag_map=mesh->TagMap;
2052 while (tag_map) {
2053 temp << tag_map->name;
2054 tag_map=tag_map->next;
2055 if (tag_map) temp << ", ";
2056 }
2057 return temp.str();
2058 }
2059
2060 } // end of namespace

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26