AcouSTO  version 2.0

◆ mshoutb()

void mshoutb ( int  icounter,
COMPLEX  cfreq 
)

Routine mshoutb(icounter,CFREQ). Produces the file 'runname.surf.fHz.msh', which can bee used to visualize the solution on the boundary. Arguments: 'icounter' is the frequency index (equal to ifreq); 'cfreq' is the complex frequency (the Laplace variable). 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
39  {
40 
41  int i,j;
42  FILE *f;
43  char fname[MAX_PATH];
44  Vec gPhinc = NULL;
45  COMPLEX val;
46 
47  int locRows;
48  int gRows;
49  int mm,nn;
50  int descl[9];
51  int descg[9];
52  int izero=0;
53  int ione=1;
54  int info,iel;
55 
56  if (0 == rank) gPhinc = calloc(modgeominfo->ncntr,sizeof(COMPLEX));
57 
58  gRows = modgeominfo->ncntr;
59  locRows = numroc_(&gRows ,&runinfo->row_block_size,&runinfo->myrow,&izero,&runinfo->nprows);
60 
61  // Gathering vectors on node 0
62  mm = gRows;// * runinfo->nprows;
63  nn = runinfo->npcols;
64  descinit_(descl , &gRows ,&ione, &runinfo->row_block_size ,&ione,&izero,&izero,&runinfo->ctxt,&locRows ,&info);
65  descinit_(descg , &mm ,&ione, &gRows ,&ione,&izero,&izero,&runinfo->ctxt,&gRows ,&info);
66  pzgemr2d_(&gRows,&ione,solution->vecPhinc ,&ione, &ione,descl, gPhinc ,&ione,&ione,descg,&runinfo->ctxt);
67 
68  if (0 != rank) return;
69  get_filename(fname,"%s-surf-%.4fHz.msh",rundetails->title,(double)HZFREQ(cfreq));
70  f = fopen(fname,"w");
71 
72  fprintf(f,"$MeshFormat\n");
73  fprintf(f,"2.2 0 %d \n",(int)sizeof(double));
74  fprintf(f,"$EndMeshFormat\n");
75 
76  fprintf(f,"$Nodes\n");
77  fprintf(f,"%d \n",modgeominfo->nnodb);
78  for(i=0;i<modgeominfo->nnodb;i++){
79  fprintf(f,"%d %f %f %f \n",i+1,(double)geometry->nodes[i].x,(double)geometry->nodes[i].y,(double)geometry->nodes[i].z);
80  }
81  fprintf(f,"$EndNodes\n");
82 
83  fprintf(f,"$Elements\n");
84  fprintf(f,"%d \n",modgeominfo->nelmb);
85  for(i=0;i<modgeominfo->nelmb;i++){
86  fprintf(f,"%d 3 2 1 1 %d %d %d %d\n",i+1,geometry->elements[i].idx[0]+1,geometry->elements[i].idx[1]+1,geometry->elements[i].idx[2]+1,geometry->elements[i].idx[3]+1);
87  }
88  fprintf(f,"$EndElements\n");
89 
90  fprintf(f,"$ElementData\n");
91  fprintf(f,"%d \n",1);
92  fprintf(f,"\"Re(Incident_Field)\"\n");
93  fprintf(f,"%d \n",1);
94  fprintf(f,"%f \n",HZFREQ(cfreq));
95  fprintf(f,"%d \n",3);
96  fprintf(f,"%d \n",icounter);
97  fprintf(f,"%d \n",1);
98  fprintf(f,"%d \n",modgeominfo->nelmb);
99  for(i=0;i<runinfo->nsymm;i++){
100  for(j=0;j<modgeominfo->ncntr;j++){
101  iel = j + i*modgeominfo->ncntr;
102  val = gPhinc[j];
103  fprintf(f,"%d %f\n",iel+1,CREAL(val));
104  }
105  }
106  fprintf(f,"$EndElementData\n");
107 
108  fprintf(f,"$ElementData\n");
109  fprintf(f,"%d \n",1);
110  fprintf(f,"\"Im(Incident_Field)\"\n");
111  fprintf(f,"%d \n",1);
112  fprintf(f,"%f \n",HZFREQ(cfreq));
113  fprintf(f,"%d \n",3);
114  fprintf(f,"%d \n",icounter);
115  fprintf(f,"%d \n",1);
116  fprintf(f,"%d \n",modgeominfo->nelmb);
117  for(i=0;i<runinfo->nsymm;i++){
118  for(j=0;j<modgeominfo->ncntr;j++){
119  iel = j + i*modgeominfo->ncntr;
120  val = gPhinc[j];
121  fprintf(f,"%d %f\n",iel+1,CIMAG(val));
122  }
123  }
124  fprintf(f,"$EndElementData\n");
125 
126  fprintf(f,"$ElementData\n");
127  fprintf(f,"%d \n",1);
128  fprintf(f,"\"Abs(Incident_Field)\"\n");
129  fprintf(f,"%d \n",1);
130  fprintf(f,"%f \n",HZFREQ(cfreq));
131  fprintf(f,"%d \n",3);
132  fprintf(f,"%d \n",icounter);
133  fprintf(f,"%d \n",1);
134  fprintf(f,"%d \n",modgeominfo->nelmb);
135  for(i=0;i<runinfo->nsymm;i++){
136  for(j=0;j<modgeominfo->ncntr;j++){
137  iel = j + i*modgeominfo->ncntr;
138  val = gPhinc[j];
139  fprintf(f,"%d %f\n",iel+1,CABS(val));
140  }
141  }
142  fprintf(f,"$EndElementData\n");
143 
144  fprintf(f,"$ElementData\n");
145  fprintf(f,"%d \n",1);
146  fprintf(f,"\"Re(Scattered_Field)\"\n");
147  fprintf(f,"%d \n",1);
148  fprintf(f,"%f \n",HZFREQ(cfreq));
149  fprintf(f,"%d \n",3);
150  fprintf(f,"%d \n",icounter);
151  fprintf(f,"%d \n",1);
152  fprintf(f,"%d \n",modgeominfo->nelmb);
153  for(i=0;i<runinfo->nsymm;i++){
154  for(j=0;j<modgeominfo->ncntr;j++){
155  iel = j + i*modgeominfo->ncntr;
156  val = solution->sol[j];
157  if(1 == modsolinfo->knw)val -= gPhinc[j];
158  fprintf(f,"%d %f\n",iel+1,CREAL(val));
159  }
160  }
161  fprintf(f,"$EndElementData\n");
162 
163  fprintf(f,"$ElementData\n");
164  fprintf(f,"%d \n",1);
165  fprintf(f,"\"Im(Scattered_Field)\"\n");
166  fprintf(f,"%d \n",1);
167  fprintf(f,"%f \n",HZFREQ(cfreq));
168  fprintf(f,"%d \n",3);
169  fprintf(f,"%d \n",icounter);
170  fprintf(f,"%d \n",1);
171  fprintf(f,"%d \n",modgeominfo->nelmb);
172  for(i=0;i<runinfo->nsymm;i++){
173  for(j=0;j<modgeominfo->ncntr;j++){
174  iel = j + i*modgeominfo->ncntr;
175  val = solution->sol[j];
176  if(1 == modsolinfo->knw)val -= gPhinc[j];
177  fprintf(f,"%d %f\n",iel+1,CIMAG(val));
178  }
179  }
180  fprintf(f,"$EndElementData\n");
181 
182  fprintf(f,"$ElementData\n");
183  fprintf(f,"%d \n",1);
184  fprintf(f,"\"Abs(Scattered_Field)\"\n");
185  fprintf(f,"%d \n",1);
186  fprintf(f,"%f \n",HZFREQ(cfreq));
187  fprintf(f,"%d \n",3);
188  fprintf(f,"%d \n",icounter);
189  fprintf(f,"%d \n",1);
190  fprintf(f,"%d \n",modgeominfo->nelmb);
191  for(i=0;i<runinfo->nsymm;i++){
192  for(j=0;j<modgeominfo->ncntr;j++){
193  iel = j + i*modgeominfo->ncntr;
194  val = solution->sol[j];
195  if(1 == modsolinfo->knw)val -= gPhinc[j];
196  fprintf(f,"%d %f\n",iel+1,CABS(val));
197  }
198  }
199  fprintf(f,"$EndElementData\n");
200 
201  fprintf(f,"$ElementData\n");
202  fprintf(f,"%d \n",1);
203  fprintf(f,"\"Re(Total_Field)\"\n");
204  fprintf(f,"%d \n",1);
205  fprintf(f,"%f \n",HZFREQ(cfreq));
206  fprintf(f,"%d \n",3);
207  fprintf(f,"%d \n",icounter);
208  fprintf(f,"%d \n",1);
209  fprintf(f,"%d \n",modgeominfo->nelmb);
210  for(i=0;i<runinfo->nsymm;i++){
211  for(j=0;j<modgeominfo->ncntr;j++){
212  iel = j + i*modgeominfo->ncntr;
213  val = solution->sol[j];
214  if(2 == modsolinfo->knw)val += gPhinc[j];
215  fprintf(f,"%d %f\n",iel+1,CREAL(val));
216  }
217  }
218  fprintf(f,"$EndElementData\n");
219 
220  fprintf(f,"$ElementData\n");
221  fprintf(f,"%d \n",1);
222  fprintf(f,"\"Im(Total_Field)\"\n");
223  fprintf(f,"%d \n",1);
224  fprintf(f,"%f \n",HZFREQ(cfreq));
225  fprintf(f,"%d \n",3);
226  fprintf(f,"%d \n",icounter);
227  fprintf(f,"%d \n",1);
228  fprintf(f,"%d \n",modgeominfo->nelmb);
229  for(i=0;i<runinfo->nsymm;i++){
230  for(j=0;j<modgeominfo->ncntr;j++){
231  iel = j + i*modgeominfo->ncntr;
232  val = solution->sol[j];
233  if(2 == modsolinfo->knw)val += gPhinc[j];
234  fprintf(f,"%d %f\n",iel+1,CIMAG(val));
235  }
236  }
237  fprintf(f,"$EndElementData\n");
238 
239  fprintf(f,"$ElementData\n");
240  fprintf(f,"%d \n",1);
241  fprintf(f,"\"Abs(Total_Field)\"\n");
242  fprintf(f,"%d \n",1);
243  fprintf(f,"%f \n",HZFREQ(cfreq));
244  fprintf(f,"%d \n",3);
245  fprintf(f,"%d \n",icounter);
246  fprintf(f,"%d \n",1);
247  fprintf(f,"%d \n",modgeominfo->nelmb);
248  for(i=0;i<runinfo->nsymm;i++){
249  for(j=0;j<modgeominfo->ncntr;j++){
250  iel = j + i*modgeominfo->ncntr;
251  val = solution->sol[j];
252  if(2 == modsolinfo->knw)val += gPhinc[j];
253  fprintf(f,"%d %f\n",iel+1,CABS(val));
254  }
255  }
256  fprintf(f,"$EndElementData\n");
257 
258 // fprintf(f,"$InterpolationScheme\n");
259 // fprintf(f,"\"Prova\" \n");
260 // fprintf(f,"$EndInterpolationScheme\n");
261 
262  fclose(f);
263 
264  VecDestroy(gPhinc);
265 
266 }
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
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 nnodb
Definition: structs.h:293
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