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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26