AcouSTO  version 2.0

◆ vtkoutmmesh()

void vtkoutmmesh ( int  icounter,
COMPLEX  cfreq 
)

Routine vtkoutmmesh(icounter, cfreq). Produces the file 'runname.mics.fHz.vtk', which can bee used to visualize the solution at the microphones array. Arguments: 'icounter' is the frequency index (equal to ifreq); 'cfreq' is the complex frequency (the Laplace variable). Data format is 'DATASET POLYDATA'. THE SUBROUDTINE CAN BE USED ONLY IF A MESH TOPOLOGY IS PROVIDED FOR MICROPHONES. COMPATIBLE WITH GMSH MICS FORMAT. Scalar fields included are:

  • Re Im and Abs of the scattered field
  • Re Im and Abs of the total field
  • Total pressure in Db
  • Insertion Loss = 20 Log [(total field)/(incident field)]
Parameters
[in]cfreqCOMPLEX frequency
[in]icounterfrequency index
44  {
45 
46  int i,i1,i2,i3,i4;
47  FILE *f;
48 // FILE *fb;
49 // FILE *ft;
50  char fname[MAX_PATH];
51  //char fnameb[MAX_PATH];
52  //char fnamet[MAX_PATH];
53  COMPLEX val;
54  Vec gPhiscm = NULL;
55  Vec gPhincm = NULL;
56  DOUBLE iloss,prat,spl;
57 
58  int locRows;
59  int gRows;
60  int mm,nn;
61  int descl[9];
62  int descg[9];
63  int izero=0;
64  int ione=1;
65  int info, celltype, cellvrtx;
66 
67  celltype=7;
68  cellvrtx=modgeominfo->nmics;
69 
70  if (0 == rank) {
71  gPhiscm = calloc(modgeominfo->nmics,sizeof(COMPLEX));
72  gPhincm = calloc(modgeominfo->nmics,sizeof(COMPLEX));
73  }
74 
75  gRows = modgeominfo->nmics;
76  locRows = numroc_(&gRows ,&runinfo->mics_block_size,&runinfo->myrow,&izero,&runinfo->nprows);
77 
78 
79  // Gathering vectors on node 0
80  mm = gRows;// * runinfo->nprows;
81  nn = runinfo->npcols;
82  descinit_(descl, &gRows,&ione, &runinfo->mics_block_size, &ione, &izero, &izero, &runinfo->ctxt, &gRows, &info);
83  descinit_(descg, &mm ,&ione, &gRows , &ione, &izero, &izero, &runinfo->ctxt, &gRows, &info);
84 
85  pzgemr2d_(&gRows,&ione,solution->vecPhium ,&ione, &ione,descl, gPhiscm ,&ione,&ione,descg,&runinfo->ctxt);
86  pzgemr2d_(&gRows,&ione,solution->vecPhincm ,&ione, &ione,descl, gPhincm ,&ione,&ione,descg,&runinfo->ctxt);
87 
88  if (0 != rank) return;
89  get_filename(fname,"%s-mics-%.4fHz.vtk",rundetails->title,(double)HZFREQ(cfreq));
90  f = fopen(fname,"w");
91 
92  fprintf(f,"# vtk DataFile Version 2.0\n");
93  fprintf(f,"Solution on the Boundary\n");
94  fprintf(f,"ASCII\n");
95 
96  fprintf(f,"DATASET UNSTRUCTURED_GRID\n");
97  fprintf(f,"POINTS %d float \n",modgeominfo->nmics);
98  for(i=0;i<modgeominfo->nmics;i++){
99  fprintf(f," %f %f %f \n",(double)geometry->mics[i].x,(double)geometry->mics[i].y,(double)geometry->mics[i].z);
100  }
101 
102  fprintf(f,"CELLS %d %d\n", modgeominfo->nemics, 5*modgeominfo->nemics);
103  for(i=0;i<modgeominfo->nemics;i++){
104  i1 = geometry->jnodm[0 + i*4];
105  i2 = geometry->jnodm[1 + i*4];
106  i3 = geometry->jnodm[2 + i*4];
107  i4 = geometry->jnodm[3 + i*4];
108  fprintf(f," %d %d %d %d %d\n",4,i1,i2,i3,i4);
109  }
110  fprintf(f,"\nCELL_TYPES %d\n", modgeominfo->nemics);
111  for(i=0;i<modgeominfo->nemics;i++){
112  if (2 == geometry->kelmm[i] || 3 == geometry->kelmm[i]) celltype = 7;
113  if (4 == geometry->kelmm[i]) celltype = 10;
114  fprintf(f,"%d\n",celltype);
115  }
116 
117  fprintf(f,"POINT_DATA %d\n",modgeominfo->nmics);
118 
119  fprintf(f,"SCALARS Re(Inc_Solution) float\n");
120  fprintf(f,"LOOKUP_TABLE default\n");
121  for(i=0;i<modgeominfo->nmics;i++){
122  val = gPhincm[i];
123  fprintf(f,f1F,CREAL(val));
124  }
125  fprintf(f,"SCALARS Im(Inc_Solution) float\n");
126  fprintf(f,"LOOKUP_TABLE default\n");
127  for(i=0;i<modgeominfo->nmics;i++){
128  val = gPhincm[i];
129  fprintf(f,f1F,CIMAG(val));
130  }
131  fprintf(f,"SCALARS Abs(Inc_Solution) float\n");
132  fprintf(f,"LOOKUP_TABLE default\n");
133  for(i=0;i<modgeominfo->nmics;i++){
134  val = gPhincm[i];
135  fprintf(f,f1F,CABS(val));
136  }
137 
138  fprintf(f,"SCALARS Re(Scat_Solution) float\n");
139  fprintf(f,"LOOKUP_TABLE default\n");
140  for(i=0;i<modgeominfo->nmics;i++){
141  val = gPhiscm[i];
142  fprintf(f,f1F,CREAL(val));
143  }
144  fprintf(f,"SCALARS Im(Scat_Solution) float\n");
145  fprintf(f,"LOOKUP_TABLE default\n");
146  for(i=0;i<modgeominfo->nmics;i++){
147  val = gPhiscm[i];
148  fprintf(f,f1F,CIMAG(val));
149  }
150  fprintf(f,"SCALARS Abs(Scat_Solution) float\n");
151  fprintf(f,"LOOKUP_TABLE default\n");
152  for(i=0;i<modgeominfo->nmics;i++){
153  val = gPhiscm[i];
154  fprintf(f,f1F,CABS(val));
155  }
156 
157  fprintf(f,"SCALARS Re(Total_Solution) float\n");
158  fprintf(f,"LOOKUP_TABLE default\n");
159  for(i=0;i<gRows;i++){
160  val = gPhiscm[i] + gPhincm[i];
161  fprintf(f,f1F,CREAL(val));
162  }
163  fprintf(f,"SCALARS Im(Total_Solution) float\n");
164  fprintf(f,"LOOKUP_TABLE default\n");
165  for(i=0;i<gRows;i++){
166  val = gPhiscm[i] + gPhincm[i];
167  fprintf(f,f1F,CIMAG(val));
168  }
169  fprintf(f,"SCALARS Abs(Total_Solution) float\n");
170  fprintf(f,"LOOKUP_TABLE default\n");
171  for(i=0;i<gRows;i++){
172  val = gPhiscm[i] + gPhincm[i];
173  fprintf(f,f1F,CABS(val));
174  }
175  fprintf(f,"SCALARS Insertion_Loss float\n");
176  fprintf(f,"LOOKUP_TABLE default\n");
177  for(i=0;i<gRows;i++){
178  prat = CABS(gPhincm[i]+gPhiscm[i])/CABS(gPhincm[i]);
179  iloss = 20*log10(1./prat);
180  fprintf(f,f1F,iloss);
181  }
182  fprintf(f,"SCALARS SPL float\n");
183  fprintf(f,"LOOKUP_TABLE default\n");
184  for(i=0;i<gRows;i++){
185  prat = 0.70710678*CABS(gPhiscm[i]+gPhincm[i]);
186  spl = 20*log10(prat) + 94.0;
187  fprintf(f,f1F,spl);
188  }
189 
190  fclose(f);
191 
192  VecDestroy(gPhincm);
193  VecDestroy(gPhiscm);
194 
195 }
struct run_info * runinfo
Definition: globals.h:34
#define Vec
Definition: types.h:63
int npcols
Definition: structs.h:113
#define VecDestroy(v)
Definition: allocation.h:111
const char * title
Definition: structs.h:135
#define CREAL(x)
Definition: functions.h:49
#define f1F
Definition: formats.h:72
int nemics
Definition: structs.h:305
int myrow
Definition: structs.h:115
int nprows
Definition: structs.h:111
void get_filename(char *filename, const char *format,...)
Definition: utils.c:59
#define CIMAG(x)
Definition: functions.h:50
int mics_block_size
Definition: structs.h:124
void pzgemr2d_(int *m, int *n, COMPLEX *A, int *ia, int *ja, int *desca, COMPLEX *B, int *ib, int *jb, int *descb, int *ictxt)
Copies a (m x n) submatrix of A from element (ia,ja) on a submatrix of B from element (ib...
int ctxt
Definition: structs.h:109
void descinit_(int *desc, int *m, int *n, int *mb, int *nb, int *ia, int *ja, int *ictxt, int *llda, int *info)
struct solution_struct * solution
Definition: globals.h:52
#define COMPLEX
Definition: types.h:48
Vec vecPhincm
Definition: structs.h:582
#define CABS(x)
Definition: functions.h:51
int nmics
Definition: structs.h:304
#define HZFREQ(x)
Definition: functions.h:28
Vec vecPhium
Definition: structs.h:586
double DOUBLE
Definition: types.h:44
int numroc_(int *m, int *mb, int *p, int *ia, int *npr)
int rank
Definition: globals.h:79
struct run_details * rundetails
Definition: globals.h:36
Geometry info structure.
Definition: structs.h:161
#define MAX_PATH
Definition: constants.h:29
#define calloc(n, size)
Definition: allocation.h:37
struct modgeom_info * modgeominfo
Definition: globals.h:38