AcouSTO  version 2.0

◆ vtkoutmu2()

void vtkoutmu2 ( int  icounter,
COMPLEX  cfreq 
)

Routine vtkoutmu2(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 UNSTRUCTURED_GRID'. AN UNSTRUCTURED GRID IS BUILT FROM THE MICS VECTOR. IT USES ONE SINGLE POLY-VERTEX CELL 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
42  {
43 
44  int i;
45  FILE *f;
46  char fname[MAX_PATH];
47  COMPLEX val;
48  Vec gPhiscm = NULL;
49  Vec gPhincm = NULL;
50  DOUBLE iloss,prat,spl;
51 
52  int locRows;
53  int gRows;
54  int mm,nn;
55  int descl[9];
56  int descg[9];
57  int izero=0;
58  int ione=1;
59  int info, celltype, cellvrtx;
60 
61  celltype=2;
62  cellvrtx=modgeominfo->nmics;
63 
64  if (0 == rank) {
65  gPhiscm = calloc(modgeominfo->nmics,sizeof(COMPLEX));
66  gPhincm = calloc(modgeominfo->nmics,sizeof(COMPLEX));
67  }
68 
69  gRows = modgeominfo->nmics;
70  locRows = numroc_(&gRows ,&runinfo->mics_block_size,&runinfo->myrow,&izero,&runinfo->nprows);
71 
72 
73  // Gathering vectors on node 0
74  mm = gRows;// * runinfo->nprows;
75  nn = runinfo->npcols;
76  descinit_(descl, &gRows,&ione, &runinfo->mics_block_size, &ione, &izero, &izero, &runinfo->ctxt, &gRows, &info);
77  descinit_(descg, &mm ,&ione, &gRows , &ione, &izero, &izero, &runinfo->ctxt, &gRows, &info);
78 
79  pzgemr2d_(&gRows,&ione,solution->vecPhium ,&ione, &ione,descl, gPhiscm ,&ione,&ione,descg,&runinfo->ctxt);
80  pzgemr2d_(&gRows,&ione,solution->vecPhincm ,&ione, &ione,descl, gPhincm ,&ione,&ione,descg,&runinfo->ctxt);
81 
82  if (0 != rank) return;
83  get_filename(fname,"%s-mics-%.4fHz.vtk",rundetails->title,(double)HZFREQ(cfreq));
84  f = fopen(fname,"w");
85 
86  fprintf(f,"# vtk DataFile Version 2.0\n");
87  fprintf(f,"Solution at Mics\n");
88  fprintf(f,"ASCII\n");
89  fprintf(f,"DATASET UNSTRUCTURED_GRID\n");
90  fprintf(f,"POINTS %d float \n",modgeominfo->nmics);
91  for(i=0;i<modgeominfo->nmics;i++){
92  fprintf(f," %f %f %f \n",(double)geometry->mics[i].x,(double)geometry->mics[i].y,(double)geometry->mics[i].z);
93  }
94 
95  fprintf(f,"CELLS %d %d\n", 1, modgeominfo->nmics+1);
96  fprintf(f,"%d",cellvrtx);
97  for(i=0;i<modgeominfo->nmics;i++){
98  fprintf(f," %d",i);
99  }
100  fprintf(f,"\nCELL_TYPES %d\n", 1);
101  fprintf(f,"%d\n",celltype);
102 
103  fprintf(f,"POINT_DATA %d\n", gRows-0);
104 
105  fprintf(f,"SCALARS Re(Inc_Solution) float\n");
106  fprintf(f,"LOOKUP_TABLE default\n");
107  for(i=0;i<modgeominfo->nmics;i++){
108  val = gPhincm[i];
109  fprintf(f,f1F,CREAL(val));
110  }
111  fprintf(f,"SCALARS Im(Inc_Solution) float\n");
112  fprintf(f,"LOOKUP_TABLE default\n");
113  for(i=0;i<modgeominfo->nmics;i++){
114  val = gPhincm[i];
115  fprintf(f,f1F,CIMAG(val));
116  }
117  fprintf(f,"SCALARS Abs(Inc_Solution) float\n");
118  fprintf(f,"LOOKUP_TABLE default\n");
119  for(i=0;i<modgeominfo->nmics;i++){
120  val = gPhincm[i];
121  fprintf(f,f1F,CABS(val));
122  }
123 
124  fprintf(f,"SCALARS Re(Scat_Solution) float\n");
125  fprintf(f,"LOOKUP_TABLE default\n");
126  for(i=0;i<modgeominfo->nmics;i++){
127  val = gPhiscm[i];
128  fprintf(f,f1F,CREAL(val));
129  }
130  fprintf(f,"SCALARS Im(Scat_Solution) float\n");
131  fprintf(f,"LOOKUP_TABLE default\n");
132  for(i=0;i<modgeominfo->nmics;i++){
133  val = gPhiscm[i];
134  fprintf(f,f1F,CIMAG(val));
135  }
136  fprintf(f,"SCALARS Abs(Scat_Solution) float\n");
137  fprintf(f,"LOOKUP_TABLE default\n");
138  for(i=0;i<modgeominfo->nmics;i++){
139  val = gPhiscm[i];
140  fprintf(f,f1F,CABS(val));
141  }
142 
143  fprintf(f,"SCALARS Re(Total_Solution) float\n");
144  fprintf(f,"LOOKUP_TABLE default\n");
145  for(i=0;i<gRows;i++){
146  val = gPhiscm[i] + gPhincm[i];
147  fprintf(f,f1F,CREAL(val));
148  }
149  fprintf(f,"SCALARS Im(Total_Solution) float\n");
150  fprintf(f,"LOOKUP_TABLE default\n");
151  for(i=0;i<gRows;i++){
152  val = gPhiscm[i] + gPhincm[i];
153  fprintf(f,f1F,CIMAG(val));
154  }
155  fprintf(f,"SCALARS Abs(Total_Solution) float\n");
156  fprintf(f,"LOOKUP_TABLE default\n");
157  for(i=0;i<gRows;i++){
158  val = gPhiscm[i] + gPhincm[i];
159  fprintf(f,f1F,CABS(val));
160  }
161  fprintf(f,"SCALARS Insertion_Loss float\n");
162  fprintf(f,"LOOKUP_TABLE default\n");
163  for(i=0;i<gRows;i++){
164  prat = CABS(gPhincm[i]+gPhiscm[i])/CABS(gPhincm[i]);
165  iloss = 20*log10(1./prat);
166  fprintf(f,f1F,iloss);
167  }
168  fprintf(f,"SCALARS SPL float\n");
169  fprintf(f,"LOOKUP_TABLE default\n");
170  for(i=0;i<gRows;i++){
171  prat = 0.70710678*CABS(gPhiscm[i]+gPhincm[i])/0.00002;
172  spl = 20*log10(prat);
173  fprintf(f,f1F,spl);
174  }
175 
176  fclose(f);
177 
178  VecDestroy(gPhincm);
179  VecDestroy(gPhiscm);
180 
181 }
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 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