AcouSTO  version 2.0

◆ vtkoutb()

void vtkoutb ( int  icounter,
COMPLEX  cfreq 
)

Routine vtkoutb(icounter,CFREQ). Produces the file 'runname.surf.fHz.vtk', which can bee used to visualize the solution on the boundary. If dirneu==1 the file contains dphi/dn on S. If dirneu==2 the file contains phi on S. Arguments: 'icounter' is the frequency index (equal to ifreq); 'cfreq' is the complex frequency. Data format is 'DATASET POLYDATA'. Scalar fields included are:

  • Re Im and Abs of the incident potential
  • Re Im and Abs of the scattered potential
  • Re Im and Abs of the total potential
Parameters
[in]cfreqCOMPLEX frequency
[in]icounterfrequency index
42  {
43 
44  int i,j,iv;
45  FILE *f;
46  char fname[MAX_PATH];
47  Vec gPhinc = NULL;
48  Vec gPsinc = NULL;
49  COMPLEX val;
50  COMPLEX vali;
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 nvrt=4;
60  int info,ive1,ive2,ive3,ive4;
61 
62  if (0 == rank) {
63  gPhinc = calloc(modgeominfo->ncntr,sizeof(COMPLEX));
64  gPsinc = calloc(modgeominfo->ncntr,sizeof(COMPLEX));
65  }
66 
67  gRows = modgeominfo->ncntr;
68  locRows = numroc_(&gRows ,&runinfo->row_block_size,&runinfo->myrow,&izero,&runinfo->nprows);
69 
70  // Gathering vectors on node 0
71  mm = gRows;// * runinfo->nprows;
72  nn = runinfo->npcols;
73  descinit_(descl , &gRows ,&ione, &runinfo->row_block_size ,&ione,&izero,&izero,&runinfo->ctxt,&locRows ,&info);
74  descinit_(descg , &mm ,&ione, &gRows ,&ione,&izero,&izero,&runinfo->ctxt,&gRows ,&info);
75  pzgemr2d_(&gRows,&ione,solution->vecPhinc ,&ione, &ione,descl, gPhinc ,&ione,&ione,descg,&runinfo->ctxt);
76  pzgemr2d_(&gRows,&ione,solution->vecPsinc ,&ione, &ione,descl, gPsinc ,&ione,&ione,descg,&runinfo->ctxt);
77 
78  if (0 != rank) return;
79  get_filename(fname,"%s-surf-%.4fHz.vtk",rundetails->title,(double)HZFREQ(cfreq));
80  f = fopen(fname,"w");
81 
82  fprintf(f,"# vtk DataFile Version 2.0\n");
83  fprintf(f,"Solution on the Boundary\n");
84  fprintf(f,"ASCII\n");
85 
86  fprintf(f,"DATASET POLYDATA\n");
87  fprintf(f,"POINTS %d float \n",modgeominfo->nelmb*4);
88  for(i=0;i<modgeominfo->nelmb;i++){
89  fprintf(f," %f %f %f \n",(double)geometry->elements[i].x[0].x,(double)geometry->elements[i].x[0].y,(double)geometry->elements[i].x[0].z);
90  fprintf(f," %f %f %f \n",(double)geometry->elements[i].x[1].x,(double)geometry->elements[i].x[1].y,(double)geometry->elements[i].x[1].z);
91  fprintf(f," %f %f %f \n",(double)geometry->elements[i].x[2].x,(double)geometry->elements[i].x[2].y,(double)geometry->elements[i].x[2].z);
92  fprintf(f," %f %f %f \n",(double)geometry->elements[i].x[3].x,(double)geometry->elements[i].x[3].y,(double)geometry->elements[i].x[3].z);
93  }
94  fprintf(f,"POLYGONS %d %d \n",modgeominfo->nelmb,5*(modgeominfo->nelmb));
95  iv= 0;
96  for(i=0;i<modgeominfo->nelmb;i++){
97  ive1 = iv++;
98  ive2 = iv++;
99  ive3 = iv++;
100  ive4 = iv++;
101  fprintf(f,"%d %d %d %d %d\n",nvrt,ive1,ive2,ive3,ive4);
102  }
103  fprintf(f,"CELL_DATA %d\n", modgeominfo->nelmb);
104  fprintf(f,"SCALARS Re(Incident_Field) float\n");
105  fprintf(f,"LOOKUP_TABLE default\n");
106  for(i=0;i<runinfo->nsymm;i++){
107  for(j=0;j<modgeominfo->ncntr;j++){
108  (1 == modsolinfo->dirneu) ? (val = gPsinc[j]) : (val = gPhinc[j]);
109  fprintf(f,f1F,CREAL(val));
110  }
111  }
112  fprintf(f,"SCALARS Im(Incident_Field) float\n");
113  fprintf(f,"LOOKUP_TABLE default\n");
114  for(i=0;i<runinfo->nsymm;i++){
115  for(j=0;j<modgeominfo->ncntr;j++){
116  (1 == modsolinfo->dirneu) ? (val = gPsinc[j]) : (val = gPhinc[j]);
117  fprintf(f,f1F,CIMAG(val));
118  }
119  }
120  fprintf(f,"SCALARS Abs(Incident_Field) float\n");
121  fprintf(f,"LOOKUP_TABLE default\n");
122  for(i=0;i<runinfo->nsymm;i++){
123  for(j=0;j<modgeominfo->ncntr;j++){
124  (1 == modsolinfo->dirneu) ? (val = gPsinc[j]) : (val = gPhinc[j]);
125  fprintf(f,f1F,CABS(val));
126  }
127  }
128  fprintf(f,"SCALARS Re(Scattered_Field) float\n");
129  fprintf(f,"LOOKUP_TABLE default\n");
130  for(i=0;i<runinfo->nsymm;i++){
131  for(j=0;j<modgeominfo->ncntr;j++){
132  val = solution->sol[j];
133  (1 == modsolinfo->dirneu) ? (vali = gPsinc[j]) : (vali = gPhinc[j]);
134  if(1 == modsolinfo->knw)val -= vali;
135  fprintf(f,f1F,CREAL(val));
136  }
137  }
138  fprintf(f,"SCALARS Im(Scattered_Field) float\n");
139  fprintf(f,"LOOKUP_TABLE default\n");
140  for(i=0;i<runinfo->nsymm;i++){
141  for(j=0;j<modgeominfo->ncntr;j++){
142  val = solution->sol[j];
143  (1 == modsolinfo->dirneu) ? (vali = gPsinc[j]) : (vali = gPhinc[j]);
144  if(1 == modsolinfo->knw)val -= vali;
145  fprintf(f,f1F,CIMAG(val));
146  }
147  }
148  fprintf(f,"SCALARS Abs(Scattered_Field) float\n");
149  fprintf(f,"LOOKUP_TABLE default\n");
150  for(i=0;i<runinfo->nsymm;i++){
151  for(j=0;j<modgeominfo->ncntr;j++){
152  val = solution->sol[j];
153  (1 == modsolinfo->dirneu) ? (vali = gPsinc[j]) : (vali = gPhinc[j]);
154  if(1 == modsolinfo->knw)val -= vali;
155  fprintf(f,f1F,CABS(val));
156  }
157  }
158  fprintf(f,"SCALARS Re(Total_Field) float\n");
159  fprintf(f,"LOOKUP_TABLE default\n");
160  for(i=0;i<runinfo->nsymm;i++){
161  for(j=0;j<modgeominfo->ncntr;j++){
162  val = solution->sol[j];
163  (1 == modsolinfo->dirneu) ? (vali = gPsinc[j]) : (vali = gPhinc[j]);
164  if(2 == modsolinfo->knw)val += vali;
165  fprintf(f,f1F,CREAL(val));
166  }
167  }
168  fprintf(f,"SCALARS Im(Total_Field) float\n");
169  fprintf(f,"LOOKUP_TABLE default\n");
170  for(i=0;i<runinfo->nsymm;i++){
171  for(j=0;j<modgeominfo->ncntr;j++){
172  val = solution->sol[j];
173  (1 == modsolinfo->dirneu) ? (vali = gPsinc[j]) : (vali = gPhinc[j]);
174  if(2 == modsolinfo->knw)val += vali;
175  fprintf(f,f1F,CIMAG(val));
176  }
177  }
178  fprintf(f,"SCALARS Abs(Total_Field) float\n");
179  fprintf(f,"LOOKUP_TABLE default\n");
180  for(i=0;i<runinfo->nsymm;i++){
181  for(j=0;j<modgeominfo->ncntr;j++){
182  val = solution->sol[j];
183  (1 == modsolinfo->dirneu) ? (vali = gPsinc[j]) : (vali = gPhinc[j]);
184  if(2 == modsolinfo->knw)val += vali;
185  fprintf(f,f1F,CABS(val));
186  }
187  }
188 
189  fclose(f);
190 
191  VecDestroy(gPhinc);
192  VecDestroy(gPsinc);
193 
194 }
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
int ncntr
Definition: structs.h:299
#define CREAL(x)
Definition: functions.h:49
#define f1F
Definition: formats.h:72
Vec vecPsinc
Definition: structs.h:580
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
COMPLEX * sol
Definition: structs.h:607
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 nelmb
Definition: structs.h:295
Vec vecPhinc
Definition: structs.h:578
int ctxt
Definition: structs.h:109
int knw
Definition: structs.h:446
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
#define CABS(x)
Definition: functions.h:51
struct modsol_info * modsolinfo
Definition: globals.h:44
#define HZFREQ(x)
Definition: functions.h:28
int nsymm
Definition: structs.h:97
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
int row_block_size
Definition: structs.h:119
#define calloc(n, size)
Definition: allocation.h:37
struct modgeom_info * modgeominfo
Definition: globals.h:38
int dirneu
Definition: structs.h:444