Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ex6DataTransfer.C
Go to the documentation of this file.
1 /* *******************************************************************
2  * Rocstar Simulation Suite *
3  * Copyright@2015, Illinois Rocstar LLC. All rights reserved. *
4  * *
5  * Illinois Rocstar LLC *
6  * Champaign, IL *
7  * www.illinoisrocstar.com *
8  * sales@illinoisrocstar.com *
9  * *
10  * License: See LICENSE file in top level of distribution package or *
11  * http://opensource.org/licenses/NCSA *
12  *********************************************************************/
13 /* *******************************************************************
14  * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, *
15  * EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES *
16  * OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND *
17  * NONINFRINGEMENT. IN NO EVENT SHALL THE CONTRIBUTORS OR *
18  * COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER *
19  * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, *
20  * ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE *
21  * USE OR OTHER DEALINGS WITH THE SOFTWARE. *
22  *********************************************************************/
23 // $Id: ex1.C,v 1.20 2008/12/06 08:43:29 mtcampbe Exp $
24 
25 #include "roccom.h"
26 #include <iostream>
27 #include <vector>
28 #include <cstdio>
29 #include <cassert>
30 #include <sstream>
31 #include <cmath>
32 
33 COM_EXTERN_MODULE( Rocface);
35 
36 void makeWindow(std::string name){
37 
38  COM_new_window(name.c_str());
39  std::string solnName = name + ".soln";
40  COM_new_attribute(solnName.c_str(), 'n', COM_DOUBLE, 3, "m/s");
41 
42  std::string compName = name + ".comp";
43  COM_new_attribute(compName.c_str(), 'n', COM_DOUBLE, 3, "m/s");
44 
45 }
46 
47 void initUnstructuredMesh( std::vector<double> &coors,
48  std::vector<int> &elmts,
49  int nrow, int ncol,
50  int rank, int nproc, int type,
51  int &nnodes,int &nelem) {
52 
53  // consider the processors as a 2*(nproc/2) grid
54  int proc_col=nproc;
55  if (nproc%2 ==0) {
56  proc_col=nproc/2;
57  }
58  else {
59  proc_col=nproc;
60  }
61 
62  int row=rank/proc_col, col=rank%proc_col;
63 
64  const double width=100., length=100.;
65 
66  std::cout << "type = " << type << std::endl;
67 
68  if(type == 0){
69  nnodes = (nrow*ncol + (nrow-1)*(ncol-1));
70  nelem = 4*(nrow-1)*(ncol-1);
71  elmts.resize(nelem*3);
72  } else if(type == 1){
73  nnodes = nrow*ncol;
74  nelem = 2*(nrow-1)*(ncol-1);
75  elmts.resize(nelem*3);
76  } else {
77  nnodes = nrow*ncol;
78  nelem = (nrow-1)*(ncol-1);
79  elmts.resize(nelem*4);
80  }
81  coors.resize(nnodes*3);
82 
83  for (int i=0; i<nrow; ++i) {
84  for (int j=0; j<ncol; ++j) {
85  coors[3*(i*ncol+j) + 0]=col*length+length/(ncol-1)*j;
86  coors[3*(i*ncol+j) + 1]=row*width+width/(nrow-1)*i;
87  coors[3*(i*ncol+j) + 2]=0;
88  }
89  }
90 
91  //generating the nodal coordinates
92  double xSpacing = 0;
93  if(ncol > 1) xSpacing = width/(ncol-1);
94  double ySpacing = 0;
95  if(nrow > 1) ySpacing = length/(nrow-1);
96 
97  if (type == 0){
98  for(int i=0; i < nrow-1; ++i){
99  for(int j=0; j < ncol-1; ++j){
100  int count = nrow*ncol;
101  int elmidx = 4*(i*(ncol-1)+j);
102 
103  coors[3*(count+i*(ncol-1)+j) + 0] = xSpacing*j + xSpacing/2 + col*length;
104  coors[3*(count+i*(ncol-1)+j) + 1] = ySpacing*i + ySpacing/2 + row*width;
105  coors[3*(count+i*(ncol-1)+j) + 2] = 0;
106 
107  // Nodes per triangle element - 'bottom' triangle
108  elmts[3*(elmidx) + 0]=i*ncol+j+1;
109  elmts[3*(elmidx) + 1]=i*(ncol-1)+j+count+1;
110  elmts[3*(elmidx) + 2]=i*ncol+j+2;
111 
112  // Nodes per triangle element - 'top' triangle
113  elmts[3*(elmidx+1) + 0]=i*ncol+j+ncol+1;
114  elmts[3*(elmidx+1) + 1]=i*ncol+j+ncol+2;
115  elmts[3*(elmidx+1) + 2]=i*(ncol-1)+j+count+1;
116 
117  // Nodes per triangle element - 'left' triangle
118  elmts[3*(elmidx+2) + 0]=i*ncol+j+1;
119  elmts[3*(elmidx+2) + 1]=i*ncol+j+ncol+1;
120  elmts[3*(elmidx+2) + 2]=i*(ncol-1)+j+count+1;
121 
122  // Nodes per triangle element - 'right' triangle
123  elmts[3*(elmidx+3) + 0]=i*ncol+j+2;
124  elmts[3*(elmidx+3) + 1]=i*(ncol-1)+j+count+1;
125  elmts[3*(elmidx+3) + 2]=i*ncol+j+ncol+2;
126 
127  }
128  }
129  }
130 
131  if (type == 1){
132  for(int i=0; i < nrow-1; ++i){
133  for(int j=0; j < ncol-1; ++j){
134  int elmidx = 2*(i*(ncol-1)+j);
135 
136  // Nodes per triangle element - 'bottom' triangle
137  elmts[3*(elmidx) + 0]=i*ncol+j+1;
138  elmts[3*(elmidx) + 1]=i*ncol+j+ncol+1;
139  elmts[3*(elmidx) + 2]=i*ncol+j+2;
140 
141  // Nodes per triangle element - 'top' triangle
142  elmts[3*(elmidx+1) + 0]=i*ncol+j+ncol+1;
143  elmts[3*(elmidx+1) + 1]=i*ncol+j+ncol+2;
144  elmts[3*(elmidx+1) + 2]=i*ncol+j+2;
145 
146  }
147  }
148  }
149 
150  if (type == 2){
151  for (int i=0; i<nrow-1; ++i) {
152  for (int j=0; j<ncol-1; ++j) {
153  int elmidx = i*(ncol-1)+j;
154 
155  elmts[4*elmidx + 0]=i*ncol+j+1;
156  elmts[4*elmidx + 1]=i*ncol+j+ncol+1;
157  elmts[4*elmidx + 2]=i*ncol+j+ncol+2;
158  elmts[4*elmidx + 3]=i*ncol+j+2;
159  }
160  }
161  }
162 }
163 
164 int main(int argc, char *argv[]){
165 
166  int comm_rank=0;
167  int comm_size=1;
168  int npanes=4;
169  int nwindows=3;
170 
171  std::vector<std::string> windowNames(nwindows);
172  std::vector< std::vector< std::vector<double> > > coords(nwindows);
173  std::vector< std::vector< std::vector<int> > > elements(nwindows);
174  std::vector< std::vector< std::vector<double> > > soln(nwindows);
175  std::vector< std::vector< std::vector<double> > > comp(nwindows);
176 
177  if ( comm_rank == 0)
178  std::cout << "Initializing Roccom and Rocface" << std::endl;
179  COM_init( &argc, &argv);
180 
181  COM_set_profiling( 1);
182 
183  COM_LOAD_MODULE_STATIC_DYNAMIC( Rocface, "RFC");
184  int RFC_clear = COM_get_function_handle("RFC.clear_overlay");
185  int RFC_read = COM_get_function_handle( "RFC.read_overlay");
186  int RFC_write = COM_get_function_handle("RFC.write_overlay");
187  int RFC_overlay = COM_get_function_handle("RFC.overlay");
188  int RFC_transfer = COM_get_function_handle( "RFC.least_squares_transfer");
189  int RFC_interp = COM_get_function_handle( "RFC.interpolate");
190  int RFC_extrap = COM_get_function_handle( "RFC.extrapolate");
191 
192 
193 
194 
195 
196  for(int i=0; i < nwindows; i++)
197  {
198  std::stringstream ss;
199  ss << "Window" << i;
200  ss >> windowNames[i];
201  if ( comm_rank == 0)
202  std::cout << "Creating window " << ss.str() << std::endl;
203  makeWindow(windowNames[i]);
204  }
205 
206  for(int i=0; i < nwindows; i++){
207  std::string &windowName(windowNames[i]);
208 
209  int nrow, ncol;
210  nrow = i+3;
211  ncol = i+2;
212  int meshType = i%3;
213  int nnodes = 0;
214  int nelems = 0;
215 
216  std::cout << windowName << ":" << std::endl;
217  std::cout << nrow << " by " << ncol << std::endl;
218 
219  coords[i].resize(npanes);
220  elements[i].resize(npanes);
221  soln[i].resize(npanes);
222  comp[i].resize(npanes);
223 
224  for(int j=0; j < npanes; j++){
225 
226  int pane_id = j+1;
227 
228  initUnstructuredMesh(coords[i][j],elements[i][j],
229  nrow,ncol,j,npanes,meshType,
230  nnodes,nelems);
231 
232  // std::cout << "coords[" << i << "][" << j << "].size() = " << coords[i][j].size()
233  // << std::endl;
234  soln[i][j].resize( coords[i][j].size() );
235  comp[i][j].resize( coords[i][j].size() );
236 
237  for(int k=0; k < soln[i][j].size(); k++){
238  soln[i][j][k] = coords[i][j][k];
239  comp[i][j][k] = -1.0;
240  // std::cout << "soln[" << i << "][" << j << "][" << k
241  // << "] = " << coords[i][j][k] + 3.0 << std::endl;
242  }
243 
244  COM_set_size( (windowName+".nc"), pane_id, nnodes);
245  COM_set_array( (windowName+".nc"), pane_id, &coords[i][j][0]);
246 
247  std::string connectivityName;
248  connectivityName = ((meshType < 2) ? (windowName+".:t3:") :
249  (windowName+".:q4:"));
250  COM_set_size(connectivityName, pane_id, nelems);
251  COM_set_array(connectivityName, pane_id, &elements[i][j][0]);
252 
253  COM_set_array((windowName+".soln"), pane_id, &soln[i][j][0]);
254  COM_set_array((windowName+".comp"), pane_id, &comp[i][j][0]);
255 
256  COM_window_init_done(windowName.c_str());
257 
258 
259  // Checking values
260 
261  // std::cout << "coords.size() = " << coords.size() << std::endl;
262  // std::cout << "coords[i].size() = " << coords[i].size() << std::endl;
263  // std::cout << "coords[i][j].size() = " << coords[i][j].size() << std::endl;
264  // std::cout << "nnodes = " << nnodes << std::endl;
265 
266  // std::cout << "elements.size() = " << elements.size() << std::endl;
267  // std::cout << "elements[i].size() = " << elements[i].size() << std::endl;
268  // std::cout << "elements[i][j].size() = " << elements[i][j].size() << std::endl;
269  // std::cout << "nelems = " << nelems << std::endl;
270 
271  // std::cout << "pane " << j << ":" << std::endl;
272  // std::cout << "coords:" << std::endl;
273  // for(int k=0; k < nnodes; k++){
274  // std::cout << coords[i][j][3*k] << " " << coords[i][j][3*k+1] << " "
275  // << coords[i][j][3*k+2] << std::endl;
276  // }
277  // std::cout << "elements:" << std::endl;
278  // if(meshType == 0 || meshType == 1){
279  // for(int k=0; k < nelems; k++){
280  // std::cout << elements[i][j][3*k] << " " << elements[i][j][3*k+1] << " "
281  // << elements[i][j][3*k+2] << std::endl;
282  // }
283  // }
284  // else if(meshType == 2){
285  // std::cout << "nelems = " << nelems << std::endl;
286  // for(int k=0; k < nelems; k++){
287  // std::cout << elements[i][j][4*k] << " " << elements[i][j][4*k+1] << " "
288  // << elements[i][j][4*k+2] << " " << elements[i][j][4*k+3]
289  // << std::endl;
290  // }
291  // }
292  // std::cout << "soln:" << std::endl;
293  // for(int k=0; k < nnodes; k++){
294  // std::cout << soln[i][j][3*k+0] << " " << soln[i][j][3*k+1] << " "
295  // << soln[i][j][3*k+2] << std::endl;
296  // }
297  // std::cout << "comp:" << std::endl;
298  // for(int k=0; k < nnodes; k++){
299  // std::cout << comp[i][j][3*k+0] << " " << comp[i][j][3*k+1] << " "
300  // << comp[i][j][3*k+2] << std::endl;
301  // }
302  // std::cout << "line " << __LINE__ << std::endl;
303  }//j Loop over panes
304  }//i Loop over windows
305 
306  // std::cout << "line " << __LINE__ << std::endl;
307 
308  int tri0_mesh = COM_get_attribute_handle("Window0.mesh");
309  int tri1_mesh = COM_get_attribute_handle("Window1.mesh");
310  int quad_mesh = COM_get_attribute_handle("Window2.mesh");
311 
312  // Perform mesh overlay
313  const char *format = "HDF";
314 
315  if( comm_size == 1){
316 
317  std::cout << "Overlaying meshes..." << std::endl;
318 
319  // COM_call_function( RFC_overlay, &tri0_mesh, &tri0_mesh);
320  // COM_call_function( RFC_write, &tri0_mesh, &tri0_mesh, "tri0", "tri_0", format);
321  // COM_call_function( RFC_clear, "Window0", "Window0");
322 
323  COM_call_function( RFC_overlay, &tri0_mesh, &tri1_mesh);
324  COM_call_function( RFC_write, &tri0_mesh, &tri1_mesh, "tri01", "tri10", format);
325  COM_call_function( RFC_clear, "Window0", "Window1");
326 
327  COM_call_function( RFC_overlay, &tri0_mesh, &quad_mesh);
328  COM_call_function( RFC_write, &tri0_mesh, &quad_mesh, "tri02", "quad20", format);
329  COM_call_function( RFC_clear, "Window0", "Window2");
330 
331  COM_call_function( RFC_overlay, &tri1_mesh, &quad_mesh);
332  COM_call_function( RFC_write, &tri1_mesh, &quad_mesh, "tri12", "quad21", format);
333  COM_call_function( RFC_clear, "Window1", "Window2");
334 
335  }
336 
337  //COM_call_function( RFC_read, &tri0_mesh, &tri0_mesh, NULL, "tri0", "tri_0", format);
338  //COM_call_function( RFC_read, &tri0_mesh, &quad_mesh, NULL, "tri0", "quad", format);
339  //COM_call_function( RFC_read, &tri1_mesh, &quad_mesh, NULL, "tri1", "quad", format);
340 
341  int tri0_soln = COM_get_attribute_handle( "Window0.soln");
342  int tri1_soln = COM_get_attribute_handle( "Window1.soln");
343  int quad_soln = COM_get_attribute_handle( "Window2.soln");
344 
345  int tri0_comp = COM_get_attribute_handle( "Window0.comp");
346  int tri1_comp = COM_get_attribute_handle( "Window1.comp");
347  int quad_comp = COM_get_attribute_handle( "Window2.comp");
348 
349  COM_call_function( RFC_read, &tri0_mesh, &tri1_mesh, NULL, "tri01", "tri10", format);
350 
351  COM_call_function( RFC_transfer, &tri1_soln, &tri0_comp);
352  int check_id = 0;
353  std::vector<std::vector< double > >::iterator paneIt = comp[check_id].begin();
354  std::vector<std::vector< double > >::iterator paneIt2 = coords[check_id].begin();
355  bool pass = true;
356  double compTol = 1e-5;
357  while(paneIt != comp[check_id].end()){
358  std::vector<double> &paneComp(*paneIt++);
359  std::vector<double> &paneCoords(*paneIt2++);
360  std::vector<double>::iterator pcIt = paneComp.begin();
361  std::vector<double>::iterator pcIt2 = paneCoords.begin();
362  while(pcIt != paneComp.end()){
363  double solnDiff = std::fabs(*pcIt - *pcIt2);
364  if(solnDiff > compTol){
365  pass = false;
366  std::cout << *pcIt << " != " << *pcIt2
367  << " (" << solnDiff << ")" << std::endl;
368  }
369  *pcIt++ = -1;
370  pcIt2++;
371  }
372  }
373  std::cout << "Nodal transfer from Triangles 1 to Triangles 0 "
374  << (pass ? "passed" : "failed") << "." << std::endl;
375 
376  COM_call_function( RFC_transfer, &tri0_soln, &tri1_comp);
377  check_id = 1;
378  paneIt = comp[check_id].begin();
379  paneIt2 = coords[check_id].begin();
380  pass = true;
381  while(paneIt != comp[check_id].end()){
382  std::vector<double> &paneComp(*paneIt++);
383  std::vector<double> &paneCoords(*paneIt2++);
384  std::vector<double>::iterator pcIt = paneComp.begin();
385  std::vector<double>::iterator pcIt2 = paneCoords.begin();
386  while(pcIt != paneComp.end()){
387  double solnDiff = std::fabs(*pcIt - *pcIt2);
388  if(solnDiff > compTol){
389  pass = false;
390  std::cout << *pcIt << " != " << *pcIt2
391  << " (" << solnDiff << ")" << std::endl;
392  }
393  *pcIt++ = -1;
394  pcIt2++;
395  }
396  }
397  std::cout << "Nodal transfer from Triangles 0 to Triangles 1 "
398  << (pass ? "passed" : "failed") << "." << std::endl;
399 
400  COM_call_function( RFC_clear, "Window0", "Window1");
401 
402  COM_call_function( RFC_read, &tri0_mesh, &quad_mesh, NULL, "tri02", "quad20", format);
403  COM_call_function( RFC_transfer, &tri0_soln, &quad_comp);
404  check_id = 2;
405  paneIt = comp[check_id].begin();
406  paneIt2 = coords[check_id].begin();
407  pass = true;
408  while(paneIt != comp[check_id].end()){
409  std::vector<double> &paneComp(*paneIt++);
410  std::vector<double> &paneCoords(*paneIt2++);
411  std::vector<double>::iterator pcIt = paneComp.begin();
412  std::vector<double>::iterator pcIt2 = paneCoords.begin();
413  while(pcIt != paneComp.end()){
414  double solnDiff = std::fabs(*pcIt - *pcIt2);
415  if(solnDiff > compTol){
416  pass = false;
417  std::cout << *pcIt << " != " << *pcIt2
418  << " (" << solnDiff << ")" << std::endl;
419  }
420  *pcIt++ = -1;
421  pcIt2++;
422  }
423  }
424  std::cout << "Nodal transfer from Triangles 0 to Quads "
425  << (pass ? "passed" : "failed") << "." << std::endl;
426 
427  COM_call_function( RFC_transfer, &quad_soln, &tri0_comp);
428  check_id = 0;
429  paneIt = comp[check_id].begin();
430  paneIt2 = coords[check_id].begin();
431  pass = true;
432  while(paneIt != comp[check_id].end()){
433  std::vector<double> &paneComp(*paneIt++);
434  std::vector<double> &paneCoords(*paneIt2++);
435  std::vector<double>::iterator pcIt = paneComp.begin();
436  std::vector<double>::iterator pcIt2 = paneCoords.begin();
437  while(pcIt != paneComp.end()){
438  double solnDiff = std::fabs(*pcIt - *pcIt2);
439  if(solnDiff > compTol){
440  pass = false;
441  std::cout << *pcIt << " != " << *pcIt2
442  << " (" << solnDiff << ")" << std::endl;
443  }
444  *pcIt++ = -1;
445  pcIt2++;
446  }
447  }
448  std::cout << "Nodal transfer from Quads to Triangles 0 "
449  << (pass ? "passed" : "failed") << "." << std::endl;
450 
451  COM_call_function( RFC_clear, "Window0", "Window2");
452 
453  COM_call_function( RFC_read, &tri1_mesh, &quad_mesh, NULL, "tri12", "quad21", format);
454  COM_call_function( RFC_transfer, &tri1_soln, &quad_comp);
455  check_id = 2;
456  paneIt = comp[check_id].begin();
457  paneIt2 = coords[check_id].begin();
458  pass = true;
459  while(paneIt != comp[check_id].end()){
460  std::vector<double> &paneComp(*paneIt++);
461  std::vector<double> &paneCoords(*paneIt2++);
462  std::vector<double>::iterator pcIt = paneComp.begin();
463  std::vector<double>::iterator pcIt2 = paneCoords.begin();
464  while(pcIt != paneComp.end()){
465  double solnDiff = std::fabs(*pcIt - *pcIt2);
466  if(solnDiff > compTol){
467  pass = false;
468  std::cout << *pcIt << " != " << *pcIt2
469  << " (" << solnDiff << ")" << std::endl;
470  }
471  *pcIt++ = -1;
472  pcIt2++;
473  }
474  }
475  std::cout << "Nodal transfer from Triangles 1 to Quads "
476  << (pass ? "passed" : "failed") << "." << std::endl;
477 
478  COM_call_function( RFC_transfer, &quad_soln, &tri1_comp);
479  check_id = 1;
480  paneIt = comp[check_id].begin();
481  paneIt2 = coords[check_id].begin();
482  pass = true;
483  while(paneIt != comp[check_id].end()){
484  std::vector<double> &paneComp(*paneIt++);
485  std::vector<double> &paneCoords(*paneIt2++);
486  std::vector<double>::iterator pcIt = paneComp.begin();
487  std::vector<double>::iterator pcIt2 = paneCoords.begin();
488  while(pcIt != paneComp.end()){
489  double solnDiff = std::fabs(*pcIt - *pcIt2);
490  if(solnDiff > compTol){
491  pass = false;
492  std::cout << *pcIt << " != " << *pcIt2
493  << " (" << solnDiff << ")" << std::endl;
494  }
495  *pcIt++ = -1;
496  pcIt2++;
497  }
498  }
499  std::cout << "Nodal transfer from Quads to Triangles 1 "
500  << (pass ? "passed" : "failed") << "." << std::endl;
501 
502  COM_call_function( RFC_clear, "Window1", "Window2");
503 
504  for(int i=0; i < nwindows; i++){
505  COM_delete_window(windowNames[i]);
506  }
507  COM_finalize();
508 
509  return 0;
510 }
const int ncol
Definition: ex1.C:95
void initUnstructuredMesh(std::vector< double > &coors, std::vector< int > &elmts, int nrow, int ncol, int rank, int nproc, int type, int &nnodes, int &nelem)
void COM_delete_window(const char *wname)
Definition: roccom_c++.h:94
j indices k indices k
Definition: Indexing.h:6
void makeWindow(std::string name)
void COM_set_size(const char *wa_str, int pane_id, int size, int ng=0)
Set sizes of for a specific attribute.
Definition: roccom_c++.h:136
This file contains the prototypes for Roccom API.
int COM_get_attribute_handle(const char *waname)
Definition: roccom_c++.h:412
double length(Vector3D *const v, int n)
void COM_finalize()
Definition: roccom_c++.h:59
Definition: Rocout.h:81
blockLoc i
Definition: read.cpp:79
int elmts[total_npanes][num_elmts][4]
Definition: ex1.C:109
void COM_window_init_done(const char *w_str, int pane_changed=true)
Definition: roccom_c++.h:102
void COM_new_window(const char *wname, MPI_Comm c=MPI_COMM_NULL)
Definition: roccom_c++.h:86
void COM_set_profiling(int i)
Definition: roccom_c++.h:550
void COM_call_function(const int wf, int argc,...)
Definition: roccom_c.C:48
void COM_set_array(const char *wa_str, int pane_id, void *addr, int strd=0, int cap=0)
Associates an array with an attribute for a specific pane.
Definition: roccom_c++.h:156
int main(int argc, char *argv[])
Definition: blastest.C:94
j indices j
Definition: Indexing.h:6
void COM_init(int *argc, char ***argv)
Definition: roccom_c++.h:57
const int nrow
Definition: ex1.C:95
void COM_new_attribute(const char *wa_str, const char loc, const int type, int ncomp, const char *unit)
Registering an attribute type.
Definition: roccom_c++.h:118
static int rank
Definition: advectest.C:66
#define COM_LOAD_MODULE_STATIC_DYNAMIC(moduleName, windowString)
Definition: roccom_basic.h:111
int COM_get_function_handle(const char *wfname)
Definition: roccom_c++.h:428
#define COM_EXTERN_MODULE(moduleName)
Definition: roccom_basic.h:116