/[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 2150 - (show annotations)
Wed Dec 10 05:57:12 2008 UTC (10 years, 8 months ago) by caltinay
File size: 84711 byte(s)
Removed unnecessary "magic name length".

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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26