CalculiX  2.13
A Free Software Three-Dimensional Structural Finite Element Program
compfluid.c File Reference
#include <unistd.h>
#include <stdio.h>
#include <math.h>
#include <stdlib.h>
#include <pthread.h>
#include "CalculiX.h"
Include dependency graph for compfluid.c:

Go to the source code of this file.

Functions

void compfluid (double **cop, ITG *nk, ITG **ipkonfp, ITG *konf, char **lakonfp, char **sidefacep, ITG *ifreestream, ITG *nfreestream, ITG *isolidsurf, ITG *neighsolidsurf, ITG *nsolidsurf, ITG *nshcon, double *shcon, ITG *nrhcon, double *rhcon, double **voldp, ITG *ntmat_, ITG *nodeboun, ITG *ndirboun, ITG *nboun, ITG *ipompc, ITG *nodempc, ITG *nmpc, ITG *ikmpc, ITG *ilmpc, ITG *ithermal, ITG *ikboun, ITG *ilboun, ITG *iturbulent, ITG *isolver, ITG *iexpl, double *ttime, double *time, double *dtime, ITG *nodeforc, ITG *ndirforc, double *xforc, ITG *nforc, ITG *nelemload, char *sideload, double *xload, ITG *nload, double *xbody, ITG *ipobody, ITG *nbody, ITG *ielmatf, char *matname, ITG *mi, ITG *ncmat_, double *physcon, ITG *istep, ITG *iinc, ITG *ibody, double *xloadold, double *xboun, double *coefmpc, ITG *nmethod, double *xforcold, double *xforcact, ITG *iamforc, ITG *iamload, double *xbodyold, double *xbodyact, double *t1old, double *t1, double *t1act, ITG *iamt1, double *amta, ITG *namta, ITG *nam, double *ampli, double *xbounold, double *xbounact, ITG *iamboun, ITG *itg, ITG *ntg, char *amname, double *t0, ITG **nelemfacep, ITG *nface, double *cocon, ITG *ncocon, double *xloadact, double *tper, ITG *jmax, ITG *jout, char *set, ITG *nset, ITG *istartset, ITG *iendset, ITG *ialset, char *prset, char *prlab, ITG *nprint, double *trab, ITG *inotr, ITG *ntrans, char *filab, char *labmpc, double *sti, ITG *norien, double *orab, char *jobnamef, char *tieset, ITG *ntie, ITG *mcs, ITG *ics, double *cs, ITG *nkon, ITG *mpcfree, ITG *memmpc_, double *fmpc, ITG *nef, ITG **inomatp, double *qfx, ITG *neifa, ITG *neiel, ITG *ielfa, ITG *ifaext, double *vfa, double *vel, ITG *ipnei, ITG *nflnei, ITG *nfaext, char *typeboun, ITG *neij, double *tincf, ITG *nactdoh, ITG *nactdohinv, ITG *ielorienf, char *jobnamec, ITG *ifatie, ITG *nstate_, double *xstate, char *orname, ITG *nblk, ITG *ielblk, ITG *istartblk, ITG *iendblk, ITG *nblket, ITG *nblkze, ITG *kon)
 

Variables

static ITG num_cpus
 

Function Documentation

◆ compfluid()

void compfluid ( double **  cop,
ITG nk,
ITG **  ipkonfp,
ITG konf,
char **  lakonfp,
char **  sidefacep,
ITG ifreestream,
ITG nfreestream,
ITG isolidsurf,
ITG neighsolidsurf,
ITG nsolidsurf,
ITG nshcon,
double *  shcon,
ITG nrhcon,
double *  rhcon,
double **  voldp,
ITG ntmat_,
ITG nodeboun,
ITG ndirboun,
ITG nboun,
ITG ipompc,
ITG nodempc,
ITG nmpc,
ITG ikmpc,
ITG ilmpc,
ITG ithermal,
ITG ikboun,
ITG ilboun,
ITG iturbulent,
ITG isolver,
ITG iexpl,
double *  ttime,
double *  time,
double *  dtime,
ITG nodeforc,
ITG ndirforc,
double *  xforc,
ITG nforc,
ITG nelemload,
char *  sideload,
double *  xload,
ITG nload,
double *  xbody,
ITG ipobody,
ITG nbody,
ITG ielmatf,
char *  matname,
ITG mi,
ITG ncmat_,
double *  physcon,
ITG istep,
ITG iinc,
ITG ibody,
double *  xloadold,
double *  xboun,
double *  coefmpc,
ITG nmethod,
double *  xforcold,
double *  xforcact,
ITG iamforc,
ITG iamload,
double *  xbodyold,
double *  xbodyact,
double *  t1old,
double *  t1,
double *  t1act,
ITG iamt1,
double *  amta,
ITG namta,
ITG nam,
double *  ampli,
double *  xbounold,
double *  xbounact,
ITG iamboun,
ITG itg,
ITG ntg,
char *  amname,
double *  t0,
ITG **  nelemfacep,
ITG nface,
double *  cocon,
ITG ncocon,
double *  xloadact,
double *  tper,
ITG jmax,
ITG jout,
char *  set,
ITG nset,
ITG istartset,
ITG iendset,
ITG ialset,
char *  prset,
char *  prlab,
ITG nprint,
double *  trab,
ITG inotr,
ITG ntrans,
char *  filab,
char *  labmpc,
double *  sti,
ITG norien,
double *  orab,
char *  jobnamef,
char *  tieset,
ITG ntie,
ITG mcs,
ITG ics,
double *  cs,
ITG nkon,
ITG mpcfree,
ITG memmpc_,
double *  fmpc,
ITG nef,
ITG **  inomatp,
double *  qfx,
ITG neifa,
ITG neiel,
ITG ielfa,
ITG ifaext,
double *  vfa,
double *  vel,
ITG ipnei,
ITG nflnei,
ITG nfaext,
char *  typeboun,
ITG neij,
double *  tincf,
ITG nactdoh,
ITG nactdohinv,
ITG ielorienf,
char *  jobnamec,
ITG ifatie,
ITG nstate_,
double *  xstate,
char *  orname,
ITG nblk,
ITG ielblk,
ITG istartblk,
ITG iendblk,
ITG nblket,
ITG nblkze,
ITG kon 
)
70  {
71 
72  /* main computational fluid dynamics routine */
73 
74  char cflag[1],*lakonf=NULL,*sideface=NULL,fncvg[132]="";
75 
76  ITG *ipointer=NULL,*mast1=NULL,*irow=NULL,*icol=NULL,*jq=NULL,
77  nzs=20000000,kode,compressible,*ifabou=NULL,*ja=NULL,
78  nfabou,im,
79  *ipkonf=NULL,*nelemface=NULL,last=0,icyclic,*iau6=NULL,
80  *inomat=NULL,ithermalref,*integerglob=NULL,iincf,
81  iconvergence=0,i,*inum=NULL,iitf,ifreefa,*neielcp=NULL,
82  *iponofa=NULL,*inofa=NULL,is,ie,*ia=NULL,*ielpropf=NULL,
83  icent=0,isti=0,iqfx=0,nfield,ndim,iorienglob,force=0,icfd=1,
84  imach=0,ikappa=0,iit,jit,iatleastonepressurebc,iturb=0,
85  *inoel=NULL,*iponoel=NULL;
86 
87  ITG nelt,isym,itol,itmax,iunit,lrgw,*igwk=NULL,ligw,ierr,*iwork=NULL,iter,
88  nsave,lenw,leniw;
89 
90  double *umfa=NULL,reltime,*doubleglob=NULL,
91  *co=NULL,*vold=NULL,*coel=NULL,*cosa=NULL,*gradvel=NULL,*gradvfa=NULL,
92  *xxn=NULL,*xxi=NULL,*xle=NULL,*xlen=NULL,*xlet=NULL,timef,dtimef,
93  *cofa=NULL,*area=NULL,*xrlfa=NULL,reltimef,ttimef,*hcfa=NULL,*cvel=NULL,
94  *au=NULL,*ad=NULL,*b=NULL,*volume=NULL,*body=NULL,*dy=NULL,
95  *advfa=NULL,*ap=NULL,*bp=NULL,*xxj=NULL,*gradkel=NULL,*gradoel=NULL,
96  *v=NULL,*velo=NULL,*veloo=NULL,*cosb=NULL,dmin,tincfguess,
97  *hel=NULL,*hfa=NULL,*auv=NULL,*adv=NULL,*bv=NULL,*sel=NULL,*gamma=NULL,
98  *gradtfa=NULL,*gradtel=NULL,*umel=NULL,*cvfa=NULL,*gradpel=NULL,
99  *eei=NULL,*ener=NULL,*thicke=NULL,*eme=NULL,c[9],*gradkfa=NULL,
100  ptimef,*stn=NULL,*qfn=NULL,*hcel=NULL,*aua=NULL,a1,a2,a3,beta,
101  *prop=NULL,*dp=NULL,*xxni=NULL,*xxnj=NULL,*xxicn=NULL,*xturb=NULL,
102  *xmach=NULL,*xkappa=NULL,urelax,*flux=NULL,velnormo[5],velnorm[5],
103  relnormt,relnormv,relnormp,relnormmax=1.e30,*temp=NULL,*auv6=NULL,
104  *adv6=NULL,*auv3=NULL,*bv3=NULL,*vela=NULL,*velaa=NULL,
105  *gradofa=NULL,betam=0.1,*gradpfa=NULL;
106 
107  double tol,*rgwk=NULL,err,*sb=NULL,*sx=NULL,*rwork=NULL,*rf=NULL;
108 
109  FILE *f1;
110 
111  co=*cop;
112  ipkonf=*ipkonfp;lakonf=*lakonfp;
113  nelemface=*nelemfacep;sideface=*sidefacep;
114  vold=*voldp;inomat=*inomatp;
115 
116 #ifdef SGI
117  ITG token;
118 #endif
119 
120  strcpy(fncvg,jobnamec);
121  strcat(fncvg,"f.cvg");
122 
123  if((f1=fopen(fncvg,"w"))==NULL){
124 // if((f1=fopen("fluidconvergence","w"))==NULL){
125  printf("*ERROR in compfluid: cannot open cvg file for writing...");
126  exit(0);
127  }
128  fprintf(f1,"temperature velocity pressure\n\n");
129 
130  urelax=.2;
131 
132  /* relative time at the end of the mechanical increment */
133 
134  reltime=(*time)/(*tper);
135 
136  /* open frd-file for fluids */
137 
138  FORTRAN(openfilefluid,(jobnamef));
139 
140  /* variables for multithreading procedure */
141 
142  ITG sys_cpus;
143  char *env,*envloc,*envsys;
144 
145  num_cpus = 0;
146  sys_cpus=0;
147 
148  /* explicit user declaration prevails */
149 
150  envsys=getenv("NUMBER_OF_CPUS");
151  if(envsys){
152  sys_cpus=atoi(envsys);
153  if(sys_cpus<0) sys_cpus=0;
154  }
155 
156  /* automatic detection of available number of processors */
157 
158  if(sys_cpus==0){
159  sys_cpus = getSystemCPUs();
160  if(sys_cpus<1) sys_cpus=1;
161  }
162 
163  /* local declaration prevails, if strictly positive */
164 
165  envloc = getenv("CCX_NPROC_CFD");
166  if(envloc){
167  num_cpus=atoi(envloc);
168  if(num_cpus<0){
169  num_cpus=0;
170  }else if(num_cpus>sys_cpus){
171  num_cpus=sys_cpus;
172  }
173  }
174 
175  /* else global declaration, if any, applies */
176 
177  env = getenv("OMP_NUM_THREADS");
178  if(num_cpus==0){
179  if (env)
180  num_cpus = atoi(env);
181  if (num_cpus < 1) {
182  num_cpus=1;
183  }else if(num_cpus>sys_cpus){
184  num_cpus=sys_cpus;
185  }
186  }
187 
188 // next line is to be inserted in a similar way for all other paralell parts
189 
190  if(*nef<num_cpus) num_cpus=*nef;
191 
192  printf(" Using up to %" ITGFORMAT " cpu(s) for CFD.\n", num_cpus);
193 
194  pthread_t tid[num_cpus];
195 
196 
197  kode=0;
198 
199  /* *iexpl==0: structure:implicit, fluid:incompressible
200  *iexpl==1: structure:implicit, fluid:compressible
201  *iexpl==2: structure:explicit, fluid:incompressible
202  *iexpl==3: structure:explicit, fluid:compressible */
203 
204  if((*iexpl==1)||(*iexpl==3)){
205  compressible=1;
206  }else{
207  compressible=0;
208  }
209 
210  /* if initial conditions are specified for the temperature,
211  it is assumed that the temperature is an unknown */
212 
213  ithermalref=*ithermal;
214  if(*ithermal==1){
215  *ithermal=2;
216  }
217 
218  /* determining the matrix structure */
219 
220  NNEW(ipointer,ITG,*nef);
221  NNEW(mast1,ITG,nzs);
222  NNEW(irow,ITG,nzs);
223  NNEW(icol,ITG,*nef);
224  NNEW(jq,ITG,*nef+1);
225 
226  mastructf(nk,konf,ipkonf,lakonf,nef,icol,jq,&mast1,&irow,
227  isolver,ipointer,&nzs,ipnei,neiel,mi);
228 
229  SFREE(ipointer);SFREE(mast1);
230 
231  NNEW(iau6,ITG,6**nef);
232  FORTRAN(create_iau6,(nef,ipnei,neiel,jq,irow,&nzs,iau6,lakonf));
233 
234  NNEW(neielcp,ITG,*nflnei);
235  FORTRAN(fill_neiel,(nef,ipnei,neiel,neielcp));
236 
237  if(compressible!=1){
238  NNEW(ia,ITG,nzs+*nef);
239  NNEW(ja,ITG,*nef+1);
240  NNEW(aua,double,nzs+*nef);
241  FORTRAN(preconvert2slapcol,(irow,ia,jq,ja,&nzs,nef));
242  }
243 
244  /* calculation geometric data */
245 
246  NNEW(coel,double,3**nef);
247  NNEW(volume,double,*nef);
248  NNEW(cosa,double,*nflnei);
249  NNEW(cosb,double,*nflnei);
250  NNEW(xxn,double,3**nflnei);
251  NNEW(xxi,double,3**nflnei);
252  NNEW(xxj,double,3**nflnei);
253  NNEW(xxni,double,3**nflnei);
254  NNEW(xxicn,double,3**nflnei);
255  NNEW(xxnj,double,3**nflnei);
256  NNEW(xle,double,*nflnei);
257  NNEW(xlen,double,*nflnei);
258  NNEW(xlet,double,*nflnei);
259  NNEW(cofa,double,3**nface);
260  NNEW(area,double,*nface);
261  NNEW(xrlfa,double,3**nface);
262  NNEW(rf,double,3**nface);
263  if(*iturbulent>0) NNEW(dy,double,*nsolidsurf);
264 
265  FORTRAN(initialcfd,(nef,ipkonf,konf,lakonf,co,coel,cofa,nface,
266  ielfa,area,ipnei,neiel,xxn,xxi,xle,xlen,xlet,xrlfa,cosa,
267  volume,neifa,xxj,cosb,&dmin,ifatie,cs,tieset,&icyclic,c,
268  neij,physcon,isolidsurf,nsolidsurf,dy,xxni,xxnj,xxicn,
269  nflnei,iturbulent,rf));
270 
271 // SFREE(xxj);
272 
273  /* storing pointers to the boundary conditions in ielfa */
274 
275  NNEW(ifabou,ITG,7**nfaext);
276  FORTRAN(applyboun,(ifaext,nfaext,ielfa,ikboun,ilboun,
277  nboun,typeboun,nelemload,nload,sideload,isolidsurf,nsolidsurf,
278  ifabou,&nfabou,nface,nodeboun,ndirboun,ikmpc,ilmpc,labmpc,nmpc,
279  nactdohinv,&compressible,&iatleastonepressurebc,ipkonf,kon,konf,
280  nblk));
281  RENEW(ifabou,ITG,nfabou);
282 
283  /* catalogueing the nodes for output purposes (interpolation at
284  the nodes */
285 
286  NNEW(iponofa,ITG,*nk);
287  NNEW(inofa,ITG,2**nface*4);
288 
289  FORTRAN(cataloguenodes,(iponofa,inofa,&ifreefa,ielfa,ifabou,ipkonf,
290  konf,lakonf,nface,nk));
291 
292  RENEW(inofa,ITG,2*ifreefa);
293 
294  /* material properties for athermal calculations
295  = calculation for which no initial thermal conditions
296  were defined */
297 
298  NNEW(umfa,double,*nface);
299  NNEW(umel,double,*nef);
300 
301  if(*ithermal==0){
302 
303  /* athermal incompressible calculations */
304 
305  /* calculating the dynamic viscosity at the element centers */
306 
307  FORTRAN(calcumel,(nef,vel,shcon,nshcon,ielmatf,ntmat_,
308  ithermal,mi,umel));
309 
310  }
311 
312 
313  if(*ithermal!=0){
314  NNEW(hcfa,double,*nface);
315  NNEW(cvel,double,*nef);
316  NNEW(cvfa,double,*nface);
317  }
318 
319  if(*nbody>0) NNEW(body,double,4**nef);
320 
321  /* v is a auxiliary field: set to zero for the calls to
322  tempload */
323 
324  NNEW(v,double,5**nk);
325 
326  /* next section is for stationary calculations */
327 
328  if(*nmethod==1){
329 
330  /* boundary conditions at the end of the mechanical
331  increment */
332 
333  FORTRAN(tempload,(xforcold,xforc,xforcact,iamforc,nforc,
334  xloadold,xload,xloadact,iamload,nload,ibody,xbody,nbody,
335  xbodyold,xbodyact,t1old,t1,t1act,iamt1,nk,amta,
336  namta,nam,ampli,time,&reltime,ttime,dtime,ithermal,nmethod,
337  xbounold,xboun,xbounact,iamboun,nboun,
338  nodeboun,ndirboun,nodeforc,ndirforc,istep,iinc,
339  co,v,itg,ntg,amname,ikboun,ilboun,nelemload,sideload,mi,
340  ntrans,trab,inotr,vold,integerglob,doubleglob,tieset,istartset,
341  iendset,ialset,ntie,nmpc,ipompc,ikmpc,ilmpc,nodempc,coefmpc,
342  ipobody,iponoel,inoel));
343 
344  /* body forces (gravity, centrifugal and Coriolis forces */
345 
346  if(*nbody>0){
347  FORTRAN(inicalcbody,(nef,body,ipobody,ibody,xbody,coel,vel,lakonf,
348  nactdohinv,&icent));
349  }
350  }
351 
352  /* extrapolating the velocity from the elements centers to the face
353  centers, thereby taking the boundary conditions into account */
354 
355  NNEW(gradvel,double,9**nef);
356  NNEW(gradvfa,double,9**nface);
357 
358  FORTRAN(extrapol_vel,(nface,ielfa,xrlfa,vel,vfa,
359  ifabou,xbounact,ipnei,nef,&icyclic,c,ifatie,xxn,gradvel,
360  gradvfa,neifa,rf,area,volume,xle,xxi,xxj,xlet,
361  coefmpc,nmpc,labmpc,ipompc,nodempc,ifaext,nfaext,nactdoh));
362 
363  /* extrapolation of the pressure at the element centers
364  to the face centers */
365 
366  NNEW(gradpel,double,3**nef);
367  NNEW(gradpfa,double,3**nface);
368 
369  FORTRAN(extrapol_pel,(nface,ielfa,xrlfa,vel,vfa,
370  ifabou,xbounact,nef,gradpel,gradpfa,neifa,rf,area,volume,
371  xle,xxi,&icyclic,xxn,ipnei,ifatie,
372  coefmpc,nmpc,labmpc,ipompc,nodempc,ifaext,nfaext,nactdoh));
373 
374  /* extrapolation of the temperature at the element centers
375  to the face centers */
376 
377  if(*ithermal>0){
378 
379  NNEW(gradtel,double,3**nef);
380  NNEW(gradtfa,double,3**nface);
381 
382  FORTRAN(extrapol_tel,(nface,ielfa,xrlfa,vel,vfa,
383  ifabou,xbounact,nef,gradtel,gradtfa,neifa,rf,area,volume,
384  xle,xxi,&icyclic,xxn,ipnei,ifatie,xload,xlet,xxj,
385  coefmpc,nmpc,labmpc,ipompc,nodempc,ifaext,nfaext,nactdoh));
386 
387  /* calculating the heat conduction at the face centers */
388 
389  FORTRAN(calchcfa,(nface,vfa,cocon,ncocon,ielmatf,ntmat_,
390  mi,ielfa,hcfa));
391 
392  if(compressible!=1){
393 
394  /* calculating the specific heat at constant volume at the
395  face centers (secant value) */
396 
397  FORTRAN(calccvfa,(nface,vfa,shcon,nshcon,ielmatf,ntmat_,
398  mi,ielfa,cvfa,physcon));
399  }else{
400 
401  /* calculating the specific heat at constant volume at the
402  face centers (secant value) */
403 
404  FORTRAN(calccvfacomp,(nface,vfa,shcon,nshcon,ielmatf,ntmat_,
405  mi,ielfa,cvfa,physcon));
406  }
407  }
408 
409  NNEW(flux,double,6**nef);
410 
411  if(compressible!=1){
412 
413  /* calculating the density at the element centers */
414 
415  FORTRAN(calcrhoel,(nef,vel,rhcon,nrhcon,ielmatf,ntmat_,
416  ithermal,mi));
417 
418  /* calculating the density at the face centers */
419 
420  FORTRAN(calcrhofa,(nface,vfa,rhcon,nrhcon,ielmatf,ntmat_,
421  ithermal,mi,ielfa));
422 
423  }else{
424 
425  /* calculating the density at the element centers */
426 
427  FORTRAN(calcrhoelcomp,(nef,vel,shcon,ielmatf,ntmat_,
428  mi));
429 
430  /* calculating the density at the face centers */
431 
432  FORTRAN(calcrhofacomp,(nface,vfa,shcon,ielmatf,ntmat_,
433  mi,ielfa,ipnei,vel,nef,flux,gradpel,gradtel,xxj,
434  &betam,xlet));
435 
436  }
437 
438  /* calculating the initial mass flux */
439 
440  FORTRAN(calcinitialflux,(area,vfa,xxn,ipnei,nef,neifa,lakonf,flux));
441 
442  /* calculating the dynamic viscosity at the face centers */
443 
444  FORTRAN(calcumfa,(nface,vfa,shcon,nshcon,ielmatf,ntmat_,
445  ithermal,mi,ielfa,umfa));
446 
447  /* extrapolation of the turbulence variables at the element centers
448  to the face centers */
449 
450  if(*iturbulent>0){
451 
452  NNEW(gradkel,double,3**nef);
453  NNEW(gradkfa,double,3**nface);
454  NNEW(gradoel,double,3**nef);
455  NNEW(gradofa,double,3**nface);
456 
457  DMEMSET(vel,7**nef,8**nef,1.);
458 
459  FORTRAN(extrapol_kel,(nface,ielfa,xrlfa,vel,vfa,
460  ifabou,xbounact,nef,gradkel,gradkfa,neifa,rf,area,volume,
461  xle,xxi,&icyclic,xxn,ipnei,ifatie,xlet,xxj,
462  coefmpc,nmpc,labmpc,ipompc,nodempc,ifaext,nfaext,nactdoh,
463  umfa,physcon));
464 
465  FORTRAN(extrapol_oel,(nface,ielfa,xrlfa,vel,vfa,
466  ifabou,xbounact,nef,gradoel,gradofa,neifa,rf,area,volume,
467  xle,xxi,&icyclic,xxn,ipnei,ifatie,xlet,xxj,
468  coefmpc,nmpc,labmpc,ipompc,nodempc,ifaext,nfaext,nactdoh,
469  umfa,physcon,dy));
470 
471  }
472 
473  /* calculating the time increment */
474 
475  FORTRAN(calcguesstincf,(nface,&dmin,vfa,umfa,cvfa,hcfa,ithermal,&tincfguess,
476  &compressible));
477 
478  /* start of the major loop */
479 
480  NNEW(advfa,double,*nface);
481  NNEW(hfa,double,3**nface);
482 
483  NNEW(ap,double,*nface);
484  NNEW(bp,double,*nface);
485 
486  NNEW(au,double,*nflnei+*nef);
487  NNEW(ad,double,*nef);
488  NNEW(b,double,*nef);
489 
490  if(*nblk!=0){
491  NNEW(auv6,double,6**nef);
492  NNEW(adv6,double,6**nef);
493  NNEW(auv3,double,3**nef);
494  NNEW(bv3,double,3**nef);
495  NNEW(vela,double,8**nef);
496  NNEW(velaa,double,8**nef);
497  }else{
498  NNEW(auv,double,*nflnei+*nef);
499  }
500 
501  NNEW(bv,double,3**nef);
502  NNEW(hel,double,3**nef);
503  NNEW(sel,double,3**nef);
504 
505  NNEW(rwork,double,*nef);
506 
507  NNEW(inum,ITG,*nk);
508 
509  NNEW(velo,double,8**nef);
510  if((compressible==0)&&(*nblk==0)) NNEW(veloo,double,8**nef);
511 
512  /* initializing velo and veloo */
513 
514  if((compressible==0)&&(*nblk==0)) memcpy(&veloo[0],&vel[0],sizeof(double)*8**nef);
515  memcpy(&velo[0],&vel[0],sizeof(double)*8**nef);
516 
517  /* check output requests */
518 
519  if((strcmp1(&filab[1914],"MACH")==0)||
520  (strcmp1(&filab[3132],"PTF")==0)||
521  (strcmp1(&filab[3219],"TTF")==0)){
522  imach=1;
523  }
524 
525  if((strcmp1(&filab[3132],"PTF")==0)||
526  (strcmp1(&filab[3219],"TTF")==0)){
527  ikappa=1;
528  }
529 
530  if(strcmp1(&filab[2088],"TURB")==0){
531  iturb=1;
532  }
533 
534  for(i=0;i<*nprint;i++){
535  if(imach==0){
536  if((strcmp1(&prlab[6*i],"MACH")==0)||
537  (strcmp1(&prlab[6*i],"PTF")==0)||
538  (strcmp1(&prlab[6*i],"TTF")==0)){
539  imach=1;
540  }
541  }
542  if(ikappa==0){
543  if((strcmp1(&prlab[6*i],"PTF")==0)||
544  (strcmp1(&prlab[6*i],"TTF")==0)){
545  ikappa=1;
546  }
547  }
548  if(iturb==0){
549  if(strcmp1(&prlab[6*i],"TURB")==0){
550  iturb=1;
551  }
552  }
553  }
554 
555  iincf=0;
556 
557  if(*tincf<=0.) *tincf=tincfguess;
558  printf("time increment for the CFD-calculations = %e\n\n",*tincf);
559 
560  ttimef=*ttime;
561  timef=*time-*dtime;
562  dtimef=*tincf;
563 
564  if(*nblk==0){
565  a1=1.5/dtimef;
566  a2=-2./dtimef;
567  a3=0.5/dtimef;
568  }else{
569  a1=1./dtimef;
570  a2=-a1;
571  }
572 
573  NNEW(temp,double,6**nef);
574  NNEW(gamma,double,*nface);
575 
576  do{
577 
578  iincf++;
579 
580  printf("fluid increment = %d\n",iincf);
581 
582  timef+=dtimef;
583  if((*time<timef)&&(*nmethod==4)){
584  dtimef-=(timef-*time);
585  timef=*time;
586  last=1;
587  beta=dtimef/(*tincf);
588  a1=(2.+beta)/(1.+beta);
589  a2=-(1.+beta)/beta;
590  a3=1./(beta*(1.+beta));
591  }
592 
593  /* starting iterations till convergence of the fluid increment */
594 
595  iit=0;
596  for(i=0;i<5;i++){velnormo[i]=0;}
597  FORTRAN(norm,(vel,velnormo,nef));
598 
599  do{
600  iit++;
601 
602  printf(" iteration = %d\n",iit);
603 
604  /* conditions for transient calculations */
605 
606  if(*nmethod==4){
607 
608  /* boundary conditions at end of fluid increment */
609 
610  FORTRAN(tempload,(xforcold,xforc,xforcact,iamforc,nforc,
611  xloadold,xload,xloadact,iamload,nload,ibody,xbody,nbody,
612  xbodyold,xbodyact,t1old,t1,t1act,iamt1,nk,amta,
613  namta,nam,ampli,&timef,&reltimef,&ttimef,&dtimef,ithermal,nmethod,
614  xbounold,xboun,xbounact,iamboun,nboun,
615  nodeboun,ndirboun,nodeforc,ndirforc,istep,iinc,
616  co,v,itg,ntg,amname,ikboun,ilboun,nelemload,sideload,mi,
617  ntrans,trab,inotr,vold,integerglob,doubleglob,tieset,istartset,
618  iendset,ialset,ntie,nmpc,ipompc,ikmpc,ilmpc,nodempc,coefmpc,
619  ipobody,iponoel,inoel));
620 
621  /* body forces (gravity, centrifugal and Coriolis forces) */
622 
623  if(*nbody>0){
624  FORTRAN(calcbody,(nef,body,ipobody,ibody,xbody,coel,vel,lakonf,
625  nactdohinv));
626  }
627 
628  }else if(icent==1){
629 
630  /* body forces (gravity, centrifugal and Coriolis forces;
631  only if centrifugal forces are active => the ensuing
632  Coriolis forces depend on the actual velocity) */
633 
634  FORTRAN(calcbody,(nef,body,ipobody,ibody,xbody,coel,vel,lakonf,
635  nactdohinv));
636  }
637 
638  /* updating of the material properties */
639 
640  if(*ithermal>0){
641 
642 
643  if(compressible!=1){
644 
645  /* calculating material data
646  density (elements+faces)
647  heat capacity at constant volume (elements+faces)
648  dynamic viscosity (elements+faces)
649  heat conduction (faces) */
650 
651  FORTRAN(materialdata_cfd,(nef,vel,shcon,nshcon,ielmatf,
652  ntmat_,mi,cvel,vfa,cocon,ncocon,physcon,cvfa,
653  ithermal,nface,umel,umfa,ielfa,hcfa,rhcon,nrhcon));
654 
655  }else{
656 
657  /* calculating material data
658  heat capacity at constant volume (elements+faces)
659  dynamic viscosity (elements+faces)
660  heat conduction (faces) */
661 
662  FORTRAN(materialdata_cfd_comp,(nef,vel,shcon,nshcon,ielmatf,
663  ntmat_,mi,cvel,vfa,cocon,ncocon,physcon,cvfa,
664  ithermal,nface,umel,umfa,ielfa,hcfa));
665  }
666 
667  }
668 
669  if(*nblk==0){
670 
671  /* filling the lhs and rhs's for the balance of momentum
672  equations */
673 
674  DMEMSET(auv,0,*nflnei+*nef,0.);
675  DMEMSET(bv,0,3**nef,0.);
676 
677  if(compressible==0){
678 
679  /* calculate gamma (Ph.D. Thesis Jasak) */
680 
681  FORTRAN(calcgamma,(nface,ielfa,vel,gradvel,gamma,xlet,xxn,xxj,
682  ipnei,&betam,nef,flux));
683 
684  mafillvmain(nef,ipnei,neifa,neiel,vfa,xxn,area,
685  auv,&auv[*nflnei],jq,irow,&nzs,bv,vel,cosa,umfa,
686  xlet,xle,gradvfa,xxi,
687  body,volume,ielfa,lakonf,ifabou,nbody,
688  &dtimef,velo,veloo,sel,xrlfa,gamma,xxj,nactdohinv,&a1,
689  &a2,&a3,flux,&icyclic,c,ifatie,iau6,xxni,xxnj,
690  iturbulent,gradvel);
691 
692  }else{
693 
694  /* calculate gamma (Ph.D. Thesis Jasak) */
695 
696 // FORTRAN(calcgamma,(nface,ielfa,vel,gradvel,gamma,xlet,xxn,xxj,
697 // ipnei,&betam,nef,flux));
698 
699  mafillvcompmain(nef,ipnei,neifa,neiel,vfa,xxn,area,
700  auv,&auv[*nflnei],jq,irow,&nzs,bv,vel,cosa,umfa,
701  xlet,xle,gradvfa,xxi,
702  body,volume,ielfa,lakonf,ifabou,nbody,
703  &dtimef,velo,veloo,sel,xrlfa,gamma,xxj,nactdohinv,&a1,
704  &a2,&a3,flux,&icyclic,c,ifatie,iau6,xxni,xxnj);
705  }
706 
707  isym=0;
708  nelt=*nflnei+*nef;
709  lrgw=131+16**nef;
710  NNEW(rgwk,double,lrgw);
711  NNEW(igwk,ITG,20);
712  for(i=0;i<*nef;i++){rwork[i]=1./auv[*nflnei+i];}
713 
714 // if(compressible==0) memcpy(&temp[*nef],&vel[*nef],sizeof(double)*3**nef);
715  memcpy(&temp[*nef],&vel[*nef],sizeof(double)*3**nef);
716 
717  /* estimate of new solution */
718 
719  FORTRAN(predgmres,(nef,&bv[0],&vel[*nef],&nelt,neielcp,ipnei,auv,
720  &isym,&itol,&tol,&itmax,&iter,
721  &err,&ierr,&iunit,sb,sx,rgwk,&lrgw,igwk,
722  &ligw,rwork,iwork));
723  if(ierr>0){
724  printf("*WARNING in compfluid: error message from predgmres (v_x)=%d\n",ierr);
725  }
726  FORTRAN(predgmres,(nef,&bv[*nef],&vel[2**nef],&nelt,neielcp,ipnei,auv,
727  &isym,&itol,&tol,&itmax,&iter,
728  &err,&ierr,&iunit,sb,sx,rgwk,&lrgw,igwk,
729  &ligw,rwork,iwork));
730  if(ierr>0){
731  printf("*WARNING in compfluid: error message from predgmres (v_y)=%d\n",ierr);
732  }
733  FORTRAN(predgmres,(nef,&bv[2**nef],&vel[3**nef],&nelt,neielcp,ipnei,auv,
734  &isym,&itol,&tol,&itmax,&iter,
735  &err,&ierr,&iunit,sb,sx,rgwk,&lrgw,igwk,
736  &ligw,rwork,iwork));
737  if(ierr>0){
738  printf("*WARNING in compfluid: error message from predgmres (v_z)=%d\n",ierr);
739  }
740  SFREE(rgwk);SFREE(igwk);
741 
742 // if(compressible==0)for(i=*nef;i<4**nef;i++){vel[i]=0.8*vel[i]+0.2*temp[i];}
743  for(i=*nef;i<4**nef;i++){vel[i]=0.8*vel[i]+0.2*temp[i];}
744 
745  }else{
746 
747  /* BLOCK structure */
748 
749  }
750 
751  /* calculating the pressure gradient at the element
752  centers */
753 
754  if(compressible==1){
755 // jit=100;
756  jit=4;
757  }else{
758 // jit=1;
759  jit=2;
760  }
761 
762  for(iitf=0;iitf<jit;iitf++){
763 
764  memcpy(&hel[0],&sel[0],sizeof(double)*(3**nef));
765 
766  /* completing hel with the neighboring velocity contributions */
767 
768  if(*nblk==0){
769  if(icyclic==0){
770  FORTRAN(complete_hel,(nef,&vel[*nef],hel,&auv[*nflnei],auv,ipnei,neiel,&nzs));
771  }else{
772  FORTRAN(complete_hel_cyclic,(nef,&vel[*nef],hel,&auv[*nflnei],auv,jq,
773  irow,ipnei,neiel,ifatie,c,lakonf,neifa,&nzs));
774  }
775  }else{
776  if(icyclic==0){
777  FORTRAN(complete_hel_blk,(vel,hel,auv6,ipnei,neiel,nef,nactdohinv));
778  }else{
779  FORTRAN(complete_hel_cyclic_blk,(vel,hel,auv6,c,ipnei,neiel,
780  neifa,ifatie,nef));
781  }
782  }
783 
784  /* generating ad and h at the face centers (advfa and hfa) */
785 
786  if(compressible!=1){
787  FORTRAN(extrapolate_ad_h,(nface,ielfa,xrlfa,&auv[*nflnei],advfa,hel,hfa,
788  &icyclic,c,ifatie));
789  }else{
790  FORTRAN(extrapolate_ad_h_comp,(nface,ielfa,xrlfa,&auv[*nflnei],advfa,hel,
791  hfa,&icyclic,c,ifatie));
792  }
793 
794  /* calculating the lhs and rhs of the equation system to determine
795  p (balance of mass) */
796 
797  DMEMSET(b,0,*nef,0.);
798 
799  if(compressible!=1){
800 
801  /* incompressible media */
802 
803  if(iitf==0){
804 
805  /* first iteration: calculating both lhs and rhs */
806 
807  DMEMSET(ad,0,*nef,0.);
808  DMEMSET(au,0,nzs,0.);
809 
810  mafillpmain(nef,lakonf,ipnei,neifa,neiel,vfa,area,
811  advfa,xlet,cosa,volume,au,ad,jq,irow,ap,
812  ielfa,ifabou,xle,b,xxn,nef,
813  &nzs,hfa,gradpel,bp,xxi,neij,xlen,cosb,
814  &iatleastonepressurebc,iau6,xxicn);
815 
816  FORTRAN(convert2slapcol,(au,ad,jq,&nzs,nef,aua));
817 
818  }else{
819 
820  /* second, third.. iteration: calculate the rhs only */
821 
822  rhspmain(nef,lakonf,ipnei,neifa,neiel,vfa,area,
823  advfa,xlet,cosa,volume,au,ad,jq,irow,ap,ielfa,ifabou,xle,
824  b,xxn,nef,&nzs,hfa,gradpel,bp,xxi,neij,xlen,
825  &iatleastonepressurebc,xxicn);
826 
827  }
828 
829 
830  nelt=nzs+*nef;
831  isym=1;
832 
833  /* next line was changed from 10 to 3 on 22.12.2016 */
834 
835  nsave=3;
836  itol=0;
837  tol=1.e-6;
838 
839  /* next line was changed from 110 to 10 on 22.12.2016 */
840 
841  itmax=10;
842  iunit=0;
843  lenw=131+17**nef+2*nelt;
844  NNEW(rgwk,double,lenw);
845  leniw=32+4**nef+2*nelt;
846  NNEW(igwk,ITG,leniw);
847 
848  memcpy(&temp[4**nef],&vel[4**nef],sizeof(double)**nef);
849 
850  FORTRAN(dslugm,(nef,&b[0],&vel[4**nef],&nelt,ia,ja,aua,
851  &isym,&nsave,&itol,&tol,&itmax,&iter,
852  &err,&ierr,&iunit,rgwk,&lenw,igwk,&leniw));
853  SFREE(rgwk);SFREE(igwk);
854 
855  for(i=4**nef;i<5**nef;i++){vel[i]=0.3*vel[i]+0.7*temp[i];}
856 
857  /* extrapolation of the pressure at the element centers
858  to the face centers */
859 
860  }else{
861 
862  /* compressible media */
863 
864  DMEMSET(au,0,*nflnei+*nef,0.);
865 
866  /* calculate gamma (Ph.D. Thesis Jasak) */
867 
868 // FORTRAN(calcgammap,(nface,ielfa,vel,gradtel,gamma,xlet,xxn,xxj,
869 // ipnei,&betam,nef,flux));
870 
871  mafillpcompmain(nef,lakonf,ipnei,neifa,neiel,vfa,area,
872  advfa,xlet,cosa,volume,au,&au[*nflnei],jq,irow,ap,
873  ielfa,ifabou,xle,b,xxn,nef,
874  &nzs,hfa,gradpel,bp,xxi,neij,xlen,cosb,
875  ielmatf,mi,&a1,&a2,&a3,velo,veloo,&dtimef,shcon,
876  ntmat_,vel,nactdohinv,xrlfa,flux,iau6,xxicn,
877  gamma);
878 
879  isym=0;
880  nelt=*nflnei+*nef;
881  lrgw=131+16**nef;
882  NNEW(rgwk,double,lrgw);
883  NNEW(igwk,ITG,20);
884  for(i=0;i<*nef;i++){rwork[i]=1./au[*nflnei+i];}
885 
886  NNEW(dp,double,*nef);
887  FORTRAN(predgmres,(nef,&b[0],dp,&nelt,neielcp,ipnei,au,
888  &isym,&itol,&tol,&itmax,&iter,
889  &err,&ierr,&iunit,sb,sx,rgwk,&lrgw,igwk,
890  &ligw,rwork,iwork));
891 
892  for(i=0;i<*nef;i++){
893 // vel[4**nef+i]+=0.2*dp[i];
894  vel[4**nef+i]+=0.3*dp[i];
895  }
896  SFREE(dp);
897 
898  SFREE(rgwk);SFREE(igwk);
899  if(ierr>0){
900  printf("*WARNING in compfluid: error message from predgmres (p)=%d\n",ierr);
901  }
902 
903  }
904 
905  FORTRAN(extrapol_pel,(nface,ielfa,xrlfa,vel,vfa,
906  ifabou,xbounact,nef,gradpel,gradpfa,neifa,rf,area,volume,
907  xle,xxi,&icyclic,xxn,ipnei,ifatie,
908  coefmpc,nmpc,labmpc,ipompc,nodempc,ifaext,nfaext,nactdoh));
909 
910  /* correction of the velocity at the element centers due
911  to the pressure change */
912 
913  FORTRAN(correctvel,(hel,&auv[*nflnei],vfa,ipnei,area,&vel[*nef],xxn,neifa,
914  lakonf,nef,nef));
915 
916  if(compressible!=0){
917  if((iitf<jit-1)&&((relnormmax>=1.e-5)||(iitf<1))){
918 
919  FORTRAN(correctvfa,(nface,ielfa,area,vfa,ap,bp,xxn,
920  ifabou,ipnei,nef,neifa,hfa,vel,xbounact,lakonf,
921  flux));
922 
923  /* calculating the lhs and rhs of the energy equation */
924 
925  DMEMSET(au,0,*nflnei+*nef,0.);
926  DMEMSET(b,0,*nef,0.);
927 
928  /* calculate gamma (Ph.D. Thesis Jasak) */
929 
930 // FORTRAN(calcgammat,(nface,ielfa,vel,gradtel,gamma,xlet,xxn,xxj,
931 // ipnei,&betam,nef,flux));
932 
933  mafilltcompmain(nef,ipnei,neifa,neiel,vfa,xxn,area,
934  au,&au[*nflnei],jq,irow,&nzs,b,vel,umel,xlet,
935  xle,gradtfa,xxi,
936  body,volume,ielfa,lakonf,ifabou,nbody,nef,
937  &dtimef,velo,veloo,cvfa,hcfa,cvel,gradvel,
938  xload,gamma,
939  xrlfa,xxj,nactdohinv,&a1,&a2,&a3,flux,iau6,
940  xxni,xxnj);
941 
942  nelt=*nflnei+*nef;
943  isym=0;
944  lrgw=131+16**nef;
945  NNEW(rgwk,double,lrgw);
946  NNEW(igwk,ITG,20);
947  for(i=0;i<*nef;i++){rwork[i]=1./au[*nflnei+i];}
948  FORTRAN(predgmres,(nef,&b[0],&vel[0],&nelt,neielcp,ipnei,au,
949  &isym,&itol,&tol,&itmax,&iter,
950  &err,&ierr,&iunit,sb,sx,rgwk,&lrgw,igwk,
951  &ligw,rwork,iwork));
952  SFREE(rgwk);SFREE(igwk);
953  if(ierr>0){
954  printf("*WARNING in compfluid: error message from predgmres (T)=%d\n",ierr);
955  }
956 
957  FORTRAN(extrapol_tel,(nface,ielfa,xrlfa,vel,vfa,
958  ifabou,xbounact,nef,gradtel,gradtfa,neifa,rf,area,volume,
959  xle,xxi,&icyclic,xxn,ipnei,ifatie,xload,xlet,xxj,
960  coefmpc,nmpc,labmpc,ipompc,nodempc,ifaext,nfaext,nactdoh));
961 
962  /* calculating the density at the element centers */
963 
964  FORTRAN(calcrhoelcomp,(nef,vel,shcon,ielmatf,ntmat_,
965  mi));
966 
967  /* calculating the density at the face centers
968  (gamma method) */
969 
970  FORTRAN(calcrhofacomp,(nface,vfa,shcon,ielmatf,ntmat_,
971  mi,ielfa,ipnei,vel,nef,flux,gradpel,gradtel,xxj,
972  &betam,xlet));
973 
974  for(i=0;i<5;i++){velnorm[i]=0;}
975  FORTRAN(norm,(vel,velnorm,nef));
976 
977  relnormt=0.;
978  relnormv=0.;
979  relnormp=0.;
980  relnormmax=0.;
981 
982  if(*ithermal!=0){
983  if(velnorm[0]/(*nef)>1.e-10){
984  relnormt=fabs(velnorm[0]-velnormo[0])/(velnormo[0]);
985  if(relnormt>relnormmax) relnormmax=relnormt;
986  }
987  }
988  if((velnorm[1]+velnorm[2]+velnorm[3])/(*nef)>1.e-10){
989  relnormv=fabs(velnorm[1]+velnorm[2]+velnorm[3]-velnormo[1]-velnormo[2]-velnormo[3])/(velnormo[1]+velnormo[2]+velnormo[3]);
990  if(relnormv>relnormmax) relnormmax=relnormv;
991  }
992  if(velnorm[4]/(*nef)>1.e-10){
993  relnormp=fabs(velnorm[4]-velnormo[4])/(velnormo[4]);
994  if(relnormp>relnormmax) relnormmax=relnormp;
995  }
996  printf("%d %11.4e %11.4e %11.4e\n",iitf,relnormt,relnormv,relnormp);
997 
998  memcpy(velnormo,velnorm,sizeof(double)*5);
999 
1000  }
1001  else{break;}
1002  }
1003 
1004  }
1005 
1006  /* adding the velocity correction at the face centers
1007  due to the balance of mass =>
1008  the resulting mass flux is correct,
1009  the face velocity vectors are not necessarily
1010  needed for energy balance, balance of momentum and
1011  the turbulence equations
1012  */
1013 
1014  FORTRAN(correctvfa,(nface,ielfa,area,vfa,ap,bp,xxn,
1015  ifabou,ipnei,nef,neifa,hfa,vel,xbounact,lakonf,
1016  flux));
1017 
1018  if(*ithermal>0){
1019 
1020  /* calculating the lhs and rhs of the energy equation */
1021 
1022  DMEMSET(ad,0,*nef,0.);
1023  DMEMSET(au,0,*nflnei+*nef,0.);
1024  DMEMSET(b,0,*nef,0.);
1025 
1026  if(compressible==0){
1027 
1028  /* calculate gamma (Ph.D. Thesis Jasak) */
1029 
1030  FORTRAN(calcgammat,(nface,ielfa,vel,gradtel,gamma,xlet,xxn,xxj,
1031  ipnei,&betam,nef,flux));
1032 
1033  mafilltmain(nef,ipnei,neifa,neiel,vfa,xxn,area,
1034  au,&au[*nflnei],jq,irow,&nzs,b,vel,umel,xlet,xle,gradtfa,xxi,
1035  body,volume,ielfa,lakonf,ifabou,nbody,nef,
1036  &dtimef,velo,veloo,cvfa,hcfa,cvel,gradvel,xload,gamma,
1037  xrlfa,xxj,nactdohinv,&a1,&a2,&a3,flux,iau6,xxni,xxnj,
1038  iturbulent);
1039 
1040  }else{
1041 
1042  /* calculate gamma (Ph.D. Thesis Jasak) */
1043 
1044 // FORTRAN(calcgammat,(nface,ielfa,vel,gradtel,gamma,xlet,xxn,xxj,
1045 // ipnei,&betam,nef,flux));
1046 
1047  mafilltcompmain(nef,ipnei,neifa,neiel,vfa,xxn,area,
1048  au,&au[*nflnei],jq,irow,&nzs,b,vel,umel,xlet,xle,gradtfa,xxi,
1049  body,volume,ielfa,lakonf,ifabou,nbody,nef,
1050  &dtimef,velo,veloo,cvfa,hcfa,cvel,gradvel,xload,gamma,
1051  xrlfa,xxj,nactdohinv,&a1,&a2,&a3,flux,iau6,xxni,xxnj);
1052  }
1053 
1054  isym=0;
1055  nelt=*nflnei+*nef;
1056  lrgw=131+16**nef;
1057  NNEW(rgwk,double,lrgw);
1058  NNEW(igwk,ITG,20);
1059  for(i=0;i<*nef;i++){rwork[i]=1./au[*nflnei+i];}
1060  FORTRAN(predgmres,(nef,&b[0],&vel[0],&nelt,neielcp,ipnei,au,
1061  &isym,&itol,&tol,&itmax,&iter,
1062  &err,&ierr,&iunit,sb,sx,rgwk,&lrgw,igwk,
1063  &ligw,rwork,iwork));
1064  SFREE(rgwk);SFREE(igwk);
1065  if(ierr>0){
1066  printf("*WARNING in compfluid: error message from predgmres (T)=%d\n",ierr);
1067  }
1068 
1069  /* extrapolation of the temperature at the element centers
1070  to the face centers */
1071 
1072  FORTRAN(extrapol_tel,(nface,ielfa,xrlfa,vel,vfa,
1073  ifabou,xbounact,nef,gradtel,gradtfa,neifa,rf,area,volume,
1074  xle,xxi,&icyclic,xxn,ipnei,ifatie,xload,xlet,xxj,
1075  coefmpc,nmpc,labmpc,ipompc,nodempc,ifaext,nfaext,nactdoh));
1076 
1077  /* recalculating the density for compressible materials */
1078 
1079  if(compressible!=0){
1080 
1081  /* calculating the density at the element centers */
1082 
1083  FORTRAN(calcrhoelcomp,(nef,vel,shcon,ielmatf,ntmat_,
1084  mi));
1085 
1086  /* calculating the density at the face centers
1087  (gamma method) */
1088 
1089  FORTRAN(calcrhofacomp,(nface,vfa,shcon,ielmatf,ntmat_,
1090  mi,ielfa,ipnei,vel,nef,flux,gradpel,gradtel,xxj,
1091  &betam,xlet));
1092 
1093  }
1094 
1095  }
1096 
1097  if(*iturbulent>0){
1098 
1099  /* calculating the lhs and rhs of the k-equation */
1100 
1101  DMEMSET(au,0,*nflnei+*nef,0.);
1102  DMEMSET(b,0,*nef,0.);
1103 
1104  /* calculate gamma (Ph.D. Thesis Jasak) */
1105 
1106  FORTRAN(calcgammak,(nface,ielfa,vel,gradkel,gamma,xlet,xxn,xxj,
1107  ipnei,&betam,nef,flux));
1108 
1109  if(compressible==0){
1110  mafillkmain(nef,ipnei,neifa,neiel,vfa,xxn,area,
1111  au,&au[*nflnei],jq,irow,&nzs,b,vel,umfa,xlet,xle,gradkfa,xxi,
1112  body,volume,ielfa,lakonf,ifabou,nbody,nef,
1113  &dtimef,velo,veloo,cvfa,hcfa,cvel,gradvel,xload,gamma,
1114  xrlfa,xxj,nactdohinv,&a1,&a2,&a3,flux,iau6,xxni,xxnj,
1115  iturbulent);
1116  }else{
1117  mafilltcompmain(nef,ipnei,neifa,neiel,vfa,xxn,area,
1118  au,&au[*nflnei],jq,irow,&nzs,b,vel,umel,xlet,xle,gradtfa,xxi,
1119  body,volume,ielfa,lakonf,ifabou,nbody,nef,
1120  &dtimef,velo,veloo,cvfa,hcfa,cvel,gradvel,xload,gamma,
1121  xrlfa,xxj,nactdohinv,&a1,&a2,&a3,flux,iau6,xxni,xxnj);
1122  }
1123 
1124  isym=0;
1125  nelt=*nflnei+*nef;
1126  lrgw=131+16**nef;
1127  NNEW(rgwk,double,lrgw);
1128  NNEW(igwk,ITG,20);
1129  for(i=0;i<*nef;i++){rwork[i]=1./au[*nflnei+i];}
1130  memcpy(&temp[0],&vel[6**nef],sizeof(double)**nef);
1131  FORTRAN(predgmres,(nef,&b[0],&temp[0],&nelt,neielcp,ipnei,au,
1132  &isym,&itol,&tol,&itmax,&iter,
1133  &err,&ierr,&iunit,sb,sx,rgwk,&lrgw,igwk,
1134  &ligw,rwork,iwork));
1135  SFREE(rgwk);SFREE(igwk);
1136  if(ierr>0){
1137  printf("*WARNING in compfluid: error message from predgmres (k)=%d\n",ierr);
1138  }
1139 
1140  /* calculating the lhs and rhs of the omega-equation */
1141 
1142  DMEMSET(au,0,*nflnei+*nef,0.);
1143  DMEMSET(b,0,*nef,0.);
1144 
1145  /* calculate gamma (Ph.D. Thesis Jasak) */
1146 
1147  FORTRAN(calcgammao,(nface,ielfa,vel,gradoel,gamma,xlet,xxn,xxj,
1148  ipnei,&betam,nef,flux));
1149 
1150  if(compressible==0){
1151  mafillomain(nef,ipnei,neifa,neiel,vfa,xxn,area,
1152  au,&au[*nflnei],jq,irow,&nzs,b,vel,umfa,xlet,xle,gradofa,xxi,
1153  body,volume,ielfa,lakonf,ifabou,nbody,nef,
1154  &dtimef,velo,veloo,cvfa,hcfa,cvel,gradvel,xload,gamma,
1155  xrlfa,xxj,nactdohinv,&a1,&a2,&a3,flux,iau6,xxni,xxnj,
1156  iturbulent,gradkel,gradoel);
1157  }else{
1158  mafilltcompmain(nef,ipnei,neifa,neiel,vfa,xxn,area,
1159  au,&au[*nflnei],jq,irow,&nzs,b,vel,umel,xlet,xle,gradtfa,xxi,
1160  body,volume,ielfa,lakonf,ifabou,nbody,nef,
1161  &dtimef,velo,veloo,cvfa,hcfa,cvel,gradvel,xload,gamma,
1162  xrlfa,xxj,nactdohinv,&a1,&a2,&a3,flux,iau6,xxni,xxnj);
1163  }
1164 
1165  isym=0;
1166  nelt=*nflnei+*nef;
1167  lrgw=131+16**nef;
1168  NNEW(rgwk,double,lrgw);
1169  NNEW(igwk,ITG,20);
1170  for(i=0;i<*nef;i++){rwork[i]=1./au[*nflnei+i];}
1171  FORTRAN(predgmres,(nef,&b[0],&vel[7**nef],&nelt,neielcp,ipnei,au,
1172  &isym,&itol,&tol,&itmax,&iter,
1173  &err,&ierr,&iunit,sb,sx,rgwk,&lrgw,igwk,
1174  &ligw,rwork,iwork));
1175  SFREE(rgwk);SFREE(igwk);
1176  if(ierr>0){
1177  printf("*WARNING in compfluid: error message from predgmres (om)=%d\n",ierr);
1178  }
1179 
1180  /* storing the updated k-values into vel */
1181 
1182  memcpy(&vel[6**nef],&temp[0],sizeof(double)**nef);
1183 
1184  /* extrapolation of the turbulence variables at the element centers
1185  to the face centers */
1186 
1187  FORTRAN(extrapol_kel,(nface,ielfa,xrlfa,vel,vfa,
1188  ifabou,xbounact,nef,gradkel,gradkfa,neifa,rf,area,volume,
1189  xle,xxi,&icyclic,xxn,ipnei,ifatie,xlet,xxj,
1190  coefmpc,nmpc,labmpc,ipompc,nodempc,ifaext,nfaext,nactdoh,
1191  umfa,physcon));
1192 
1193  FORTRAN(extrapol_oel,(nface,ielfa,xrlfa,vel,vfa,
1194  ifabou,xbounact,nef,gradoel,gradofa,neifa,rf,area,volume,
1195  xle,xxi,&icyclic,xxn,ipnei,ifatie,xlet,xxj,
1196  coefmpc,nmpc,labmpc,ipompc,nodempc,ifaext,nfaext,nactdoh,
1197  umfa,physcon,dy));
1198 
1199  }
1200 
1201  /* extrapolating the velocity from the elements centers to the face
1202  centers, thereby taking the boundary conditions into account */
1203 
1204  FORTRAN(extrapol_vel,(nface,ielfa,xrlfa,vel,vfa,
1205  ifabou,xbounact,ipnei,nef,&icyclic,c,ifatie,xxn,gradvel,
1206  gradvfa,neifa,rf,area,volume,xle,xxi,xxj,xlet,
1207  coefmpc,nmpc,labmpc,ipompc,nodempc,ifaext,nfaext,nactdoh));
1208 
1209 // FORTRAN(writevfa,(vfa,nface,nactdohinv,ielfa));
1210 
1211  /* end subiterations */
1212 
1213  for(i=0;i<5;i++){velnorm[i]=0;}
1214  FORTRAN(norm,(vel,velnorm,nef));
1215 
1216  relnormt=0.;
1217  relnormv=0.;
1218  relnormp=0.;
1219  relnormmax=0.;
1220 
1221  if(*ithermal!=0){
1222  if(velnorm[0]/(*nef)>1.e-10){
1223  relnormt=fabs(velnorm[0]-velnormo[0])/(velnormo[0]);
1224  if(relnormt>relnormmax) relnormmax=relnormt;
1225  }
1226  }
1227  if((velnorm[1]+velnorm[2]+velnorm[3])/(*nef)>1.e-10){
1228  relnormv=fabs(velnorm[1]+velnorm[2]+velnorm[3]-velnormo[1]-velnormo[2]-velnormo[3])/(velnormo[1]+velnormo[2]+velnormo[3]);
1229  if(relnormv>relnormmax) relnormmax=relnormv;
1230  }
1231  if(velnorm[4]/(*nef)>1.e-10){
1232  relnormp=fabs(velnorm[4]-velnormo[4])/(velnormo[4]);
1233  if(relnormp>relnormmax) relnormmax=relnormp;
1234  }
1235  if(iit==1){
1236  fprintf(f1,"%11.4e %11.4e %11.4e\n",relnormt,relnormv,relnormp);
1237  }
1238 
1239  memcpy(velnormo,velnorm,sizeof(double)*5);
1240 
1241  if((*nmethod==1)&&(compressible!=1)){
1242 
1243  /* steady state incompressible flow:
1244  calculate the velocity only once in each increment */
1245 
1246  if(relnormmax<1.e-10) iconvergence=1;
1247 // if(relnormmax<1.e-5) iconvergence=1;
1248  break;
1249  }else{
1250 
1251  /* compressible flow:
1252  calculate the velocity only once in each increment */
1253 
1254  if((compressible==1)&&(iit==1))break;
1255 
1256  /* incompressible transient flow:
1257  calculate the velocity repeatedly in each increment */
1258 
1259  if(relnormmax<1.e-3)break;
1260  }
1261 
1262  }while(1);
1263 
1264  if(((iincf/jout[1])*jout[1]==iincf)||(iconvergence==1)||
1265  (iincf==jmax[1])){
1266 
1267  /* calculating the stress and the heat flow at the
1268  integration points, if requested */
1269 
1270  if((strcmp1(&filab[3306],"SF ")==0)||
1271  (strcmp1(&filab[3480],"SVF ")==0))isti=1;
1272  if(strcmp1(&filab[3393],"HFLF")==0)iqfx=1;
1273  for(i=0;i<*nprint;i++){
1274  if(strcmp1(&prlab[6*i],"SVF")==0) isti=1;
1275  if(strcmp1(&prlab[6*i],"HFLF")==0)iqfx=1;
1276  }
1277 
1278  /* calculating the heat conduction at the element centers */
1279 
1280  if(iqfx==1){
1281  NNEW(hcel,double,*nef);
1282  FORTRAN(calchcel,(vel,cocon,ncocon,ielmatf,ntmat_,mi,
1283  hcel,nef));
1284  }
1285 
1286  /* calculating the stress and/or the heat flux at the
1287  element centers */
1288 
1289  if((isti==1)||(iqfx==1)){
1290  FORTRAN(calcstressheatflux,(sti,umel,gradvel,qfx,hcel,
1291  gradtel,nef,&isti,&iqfx,mi));
1292  if(iqfx==1)SFREE(hcel);
1293  }
1294 
1295  /* extrapolating the stresses */
1296 
1297  if((strcmp1(&filab[3306],"SF ")==0)||
1298  (strcmp1(&filab[3480],"SVF ")==0)){
1299  nfield=6;
1300  ndim=6;
1301  if((*norien>0)&&
1302  ((strcmp1(&filab[3311],"L")==0)||(strcmp1(&filab[3485],"L")==0))){
1303  iorienglob=1;
1304  }else{
1305  iorienglob=0;
1306  }
1307  strcpy1(&cflag[0],&filab[2962],1);
1308  NNEW(stn,double,6**nk);
1309  FORTRAN(extrapolate,(sti,stn,ipkonf,inum,konf,lakonf,
1310  &nfield,nk,nef,mi,&ndim,orab,ielorienf,co,&iorienglob,
1311  cflag,vold,&force,ielmatf,thicke,ielpropf,prop));
1312  }
1313 
1314  /* extrapolating the heat flow */
1315 
1316 
1317  if(strcmp1(&filab[3393],"HFLF")==0){
1318  nfield=3;
1319  ndim=3;
1320  if((*norien>0)&&(strcmp1(&filab[3398],"L")==0)){
1321  iorienglob=1;
1322  }else{
1323  iorienglob=0;
1324  }
1325  strcpy1(&cflag[0],&filab[3049],1);
1326  NNEW(qfn,double,3**nk);
1327  FORTRAN(extrapolate,(qfx,qfn,ipkonf,inum,konf,lakonf,
1328  &nfield,nk,nef,mi,&ndim,orab,ielorienf,co,&iorienglob,
1329  cflag,vold,&force,ielmatf,thicke,ielpropf,prop));
1330  }
1331 
1332  /* extrapolating the facial values of the static temperature
1333  and/or the velocity and/or the static pressure to the nodes */
1334 
1335  if(imach){NNEW(xmach,double,*nk);}
1336  if(ikappa){NNEW(xkappa,double,*nk);}
1337  if(iturb){NNEW(xturb,double,2**nk);}
1338 
1339  FORTRAN(extrapolatefluid,(nk,iponofa,inofa,inum,vfa,vold,ielfa,
1340  ithermal,&imach,&ikappa,xmach,xkappa,shcon,nshcon,ntmat_,
1341  ielmatf,physcon,mi,&iturb,xturb));
1342 
1343  /* storing the results in dat-format */
1344 
1345  ptimef=ttimef+timef;
1346  FORTRAN(printoutfluid,(set,nset,istartset,iendset,ialset,nprint,
1347  prlab,prset,ipkonf,lakonf,sti,eei,
1348  xstate,ener,mi,nstate_,co,konf,qfx,
1349  &ptimef,trab,inotr,ntrans,orab,ielorienf,
1350  norien,vold,ielmatf,
1351  thicke,eme,xturb,physcon,nactdoh,
1352  ielpropf,prop,xkappa,xmach,ithermal,
1353  orname));
1354 
1355  /* thermal flux and drag: storage in dat-format */
1356 
1357  FORTRAN(printoutface,(co,rhcon,nrhcon,ntmat_,vold,shcon,nshcon,
1358  cocon,ncocon,&compressible,istartset,iendset,ipkonf,
1359  lakonf,konf,
1360  ialset,prset,&ptimef,nset,set,nprint,prlab,ielmatf,mi,
1361  ithermal,nactdoh,&icfd,time,stn));
1362 
1363  /* storing the results in frd-format */
1364 
1365  FORTRAN(frdfluid,(co,nk,konf,ipkonf,lakonf,nef,vold,&kode,&timef,ielmatf,
1366  matname,filab,inum,ntrans,inotr,trab,mi,istep,
1367  stn,qfn,nactdohinv,xmach,xkappa,physcon,xturb));
1368 
1369 // FORTRAN(writevfa,(vfa,nface,nactdohinv,ielfa));
1370 
1371  if((strcmp1(&filab[3306],"SF ")==0)||
1372  (strcmp1(&filab[3480],"SVF ")==0)){SFREE(stn);}
1373  if(strcmp1(&filab[3393],"HFLF")==0){SFREE(qfn);}
1374 
1375  if(imach){SFREE(xmach);}
1376  if(ikappa){SFREE(xkappa);}
1377  if(iturb){SFREE(xturb);}
1378 
1379  }
1380 
1381  if(iincf==jmax[1]){
1382  printf("*INFO: maximum number of fluid increments reached\n\n");
1383  fclose(f1);
1384  FORTRAN(stop,());
1385  }
1386  if(last==1){
1387  printf("*INFO: mechanical time increment reached: time=%e\n\n",*dtime);
1388  fclose(f1);
1389  FORTRAN(stop,());
1390  }
1391  if(iconvergence==1){
1392  printf("*INFO: steady state reached\n\n");
1393  fclose(f1);
1394  FORTRAN(stop,());
1395  }
1396 
1397 
1398  if((compressible==0)&&(*nblk==0)) memcpy(&veloo[0],&velo[0],sizeof(double)*8**nef);
1399  memcpy(&velo[0],&vel[0],sizeof(double)*8**nef);
1400 
1401  }while(1);
1402 
1403  FORTRAN(closefilefluid,());
1404 
1405  SFREE(flux);
1406 
1407  if(compressible!=1){SFREE(ia);SFREE(ja);SFREE(aua);}
1408 
1409  SFREE(irow);SFREE(icol);SFREE(jq);SFREE(iau6);SFREE(neielcp);
1410 
1411  SFREE(coel);SFREE(cosa);SFREE(xxn);SFREE(xxi);SFREE(xle);SFREE(xlen);
1412  SFREE(xlet);SFREE(cofa);SFREE(area);SFREE(xrlfa);SFREE(volume);
1413  SFREE(cosb);SFREE(xxni);SFREE(xxnj);SFREE(xxicn);SFREE(xxj);
1414  SFREE(rf);
1415  if(*iturbulent>0) SFREE(dy);
1416 
1417  SFREE(ifabou);SFREE(umfa);SFREE(umel);
1418 
1419  SFREE(gradvel);SFREE(gradvfa);SFREE(au);SFREE(ad);SFREE(b);SFREE(advfa);
1420  SFREE(ap);SFREE(bp);SFREE(gradpel);SFREE(rwork);
1421  SFREE(hfa);SFREE(hel);SFREE(adv);SFREE(bv);SFREE(sel);
1422  if(*nblk!=0){
1423  SFREE(auv6);SFREE(adv6);SFREE(auv3);SFREE(bv3);
1424  SFREE(vela);SFREE(velaa);
1425  }else{
1426  SFREE(auv);
1427  }
1428 
1429  if(*ithermal>0){
1430  SFREE(gradtel);SFREE(gradtfa);SFREE(hcfa);SFREE(cvel);SFREE(cvfa);
1431  }
1432 
1433  if(*iturbulent>0){
1434  SFREE(gradkel);SFREE(gradkfa);SFREE(gradoel);SFREE(gradofa);
1435  }
1436 
1437  SFREE(inum);SFREE(v);SFREE(velo);
1438  if((compressible==0)&&(*nblk==0)) SFREE(veloo);
1439 
1440  SFREE(iponofa);SFREE(inofa);
1441 
1442  if(*nbody>0) SFREE(body);
1443 
1444  *ithermal=ithermalref;
1445 
1446  SFREE(temp);SFREE(gamma);
1447 
1448  SFREE(gradpfa);
1449 
1450  return;
1451 
1452 }
#define ITGFORMAT
Definition: CalculiX.h:52
subroutine calchcel(vel, cocon, ncocon, ielmat, ntmat_, mi, hcel, nef)
Definition: calchcel.f:21
subroutine extrapolate_ad_h_comp(nface, ielfa, xrlfa, adv, advfa, hel, hfa, icyclic, c, ifatie)
Definition: extrapolate_ad_h_comp.f:21
subroutine convert2slapcol(au, ad, jq, nzs, nef, aua)
Definition: convert2slapcol.f:26
void mafilltmain(ITG *nef, ITG *ipnei, ITG *neifa, ITG *neiel, double *vfa, double *xxn, double *area, double *au, double *ad, ITG *jq, ITG *irow, ITG *nzs, double *b, double *vel, double *umel, double *xlet, double *xle, double *gradtfa, double *xxi, double *body, double *volume, ITG *ielfa, char *lakonf, ITG *ifabou, ITG *nbody, ITG *neq, double *dtimef, double *velo, double *veloo, double *cpfa, double *hcfa, double *cvel, double *gradvel, double *xload, double *gammat, double *xrlfa, double *xxj, ITG *nactdohinv, double *a1, double *a2, double *a3, double *flux, ITG *iau6, double *xxni, double *xxnj, ITG *iturbulent)
Definition: mafilltmain.c:35
subroutine calcguesstincf(nface, dmin, vfa, umfa, cvfa, hcfa, ithermal, tincfguess, compressible)
Definition: calcguesstincf.f:21
subroutine extrapolatefluid(nk, iponofa, inofa, inum, vfa, v, ielfa, ithermal, imach, ikappa, xmach, xkappa, shcon, nshcon, ntmat_, ielmat, physcon, mi, iturb, xturb)
Definition: extrapolatefluid.f:22
subroutine calccvfa(nface, vfa, shcon, nshcon, ielmat, ntmat_, mi, ielfa, cvfa, physcon)
Definition: calccvfa.f:21
subroutine complete_hel(nef, bv, hel, adv, auv, ipnei, neiel, nzs)
Definition: complete_hel.f:26
subroutine calcbody(nef, body, ipobody, ibody, xbody, coel, vel, lakon, nactdohinv)
Definition: calcbody.f:21
subroutine extrapol_oel(nface, ielfa, xrlfa, vel, vfa, ifabou, xbounact, nef, gradoel, gradofa, neifa, rf, area, volume, xle, xxi, icyclic, xxn, ipnei, ifatie, xlet, xxj, coefmpc, nmpc, labmpc, ipompc, nodempc, ifaext, nfaext, nactdoh, umfa, physcon, dy)
Definition: extrapol_oel.f:24
subroutine inicalcbody(nef, body, ipobody, ibody, xbody, coel, vel, lakon, nactdohinv, icent)
Definition: inicalcbody.f:21
subroutine extrapol_vel(nface, ielfa, xrlfa, vel, vfa, ifabou, xbounact, ipnei, nef, icyclic, c, ifatie, xxn, gradvel, gradvfa, neifa, rf, area, volume, xle, xxi, xxj, xlet, coefmpc, nmpc, labmpc, ipompc, nodempc, ifaext, nfaext, nactdoh)
Definition: extrapol_vel.f:23
subroutine printoutfluid(set, nset, istartset, iendset, ialset, nprint, prlab, prset, ipkonf, lakonf, sti, eei, xstate, ener, mi, nstate_, co, konf, qfx, ttime, trab, inotr, ntrans, orab, ielorienf, norien, vold, ielmatf, thicke, eme, xturb, physcon, nactdoh, ielpropf, prop, xkappa, xmach, ithermal, orname)
Definition: printoutfluid.f:25
subroutine cataloguenodes(iponofa, inofa, ifreefa, ielfa, ifabou, ipkon, konf, lakon, nface, nk)
Definition: cataloguenodes.f:21
subroutine calcgammat(nface, ielfa, vel, gradtel, gamma, xlet, xxn, xxj, ipnei, betam, nef, flux)
Definition: calcgammat.f:21
subroutine calcumel(nef, vel, shcon, nshcon, ielmat, ntmat_, ithermal, mi, umel)
Definition: calcumel.f:21
subroutine calcrhofa(nface, vfa, rhcon, nrhcon, ielmat, ntmat_, ithermal, mi, ielfa)
Definition: calcrhofa.f:21
void mafillpcompmain(ITG *ne, char *lakonf, ITG *ipnei, ITG *neifa, ITG *neiel, double *vfa, double *area, double *adfa, double *xlet, double *cosa, double *volume, double *au, double *ad, ITG *jq, ITG *irow, double *ap, ITG *ielfa, ITG *ifabou, double *xle, double *b, double *xxn, ITG *neq, ITG *nzs, double *hfa, double *gradpel, double *bp, double *xxi, ITG *neij, double *xlen, double *cosb, ITG *ielmatf, ITG *mi, double *a1, double *a2, double *a3, double *velo, double *veloo, double *dtimef, double *shcon, ITG *ntmat_, double *vel, ITG *nactdohinv, double *xrlfa, double *flux, ITG *iau6, double *xxicn, double *gamma)
Definition: mafillpcompmain.c:36
subroutine complete_hel_blk(vel, hel, auv6, ipnei, neiel, nef, nactdohinv)
Definition: complete_hel_blk.f:27
ITG strcpy1(char *s1, const char *s2, ITG length)
Definition: strcpy1.c:24
void FORTRAN(actideacti,(char *set, ITG *nset, ITG *istartset, ITG *iendset, ITG *ialset, char *objectset, ITG *ipkon, ITG *ibject, ITG *ne))
subroutine flux(node1, node2, nodem, nelem, lakon, kon, ipkon, nactdog, identity, ielprop, prop, kflag, v, xflow, f, nodef, idirf, df, cp, R, rho, physcon, g, co, dvi, numf, vold, set, shcon, nshcon, rhcon, nrhcon, ntmat_, mi, ider, ttime, time, iaxial)
Definition: flux.f:24
void mafillkmain(ITG *nef, ITG *ipnei, ITG *neifa, ITG *neiel, double *vfa, double *xxn, double *area, double *au, double *ad, ITG *jq, ITG *irow, ITG *nzs, double *b, double *vel, double *umel, double *xlet, double *xle, double *gradkfa, double *xxi, double *body, double *volume, ITG *ielfa, char *lakonf, ITG *ifabou, ITG *nbody, ITG *neq, double *dtimef, double *velo, double *veloo, double *cpfa, double *hcfa, double *cvel, double *gradvel, double *xload, double *gammat, double *xrlfa, double *xxj, ITG *nactdohinv, double *a1, double *a2, double *a3, double *flux, ITG *iau6, double *xxni, double *xxnj, ITG *iturbulent)
Definition: mafillkmain.c:35
ITG strcmp1(const char *s1, const char *s2)
Definition: strcmp1.c:24
subroutine calcinitialflux(area, vfa, xxn, ipnei, nef, neifa, lakonf, flux)
Definition: calcinitialflux.f:21
subroutine frdfluid(co, nk, konf, ipkonf, lakonf, nef, vold, kode, time, ielmat, matname, filab, inum, ntrans, inotr, trab, mi, istep, stn, qfn, nactdohinv, xmach, xkappa, physcon, xturb)
Definition: frdfluid.f:22
subroutine norm(vel, velnorm, nef)
Definition: norm.f:20
subroutine openfilefluid(jobname)
Definition: openfilefluid.f:20
ITG getSystemCPUs()
Definition: getSystemCPUs.c:40
subroutine tempload(xforcold, xforc, xforcact, iamforc, nforc, xloadold, xload, xloadact, iamload, nload, ibody, xbody, nbody, xbodyold, xbodyact, t1old, t1, t1act, iamt1, nk, amta, namta, nam, ampli, time, reltime, ttime, dtime, ithermal, nmethod, xbounold, xboun, xbounact, iamboun, nboun, nodeboun, ndirboun, nodeforc, ndirforc, istep, iinc, co, vold, itg, ntg, amname, ikboun, ilboun, nelemload, sideload, mi, ntrans, trab, inotr, veold, integerglob, doubleglob, tieset, istartset, iendset, ialset, ntie, nmpc, ipompc, ikmpc, ilmpc, nodempc, coefmpc, ipobody, iponoel, inoel)
Definition: tempload.f:29
#define DMEMSET(a, b, c, d)
Definition: CalculiX.h:45
subroutine extrapolate_ad_h(nface, ielfa, xrlfa, adv, advfa, hel, hfa, icyclic, c, ifatie)
Definition: extrapolate_ad_h.f:21
subroutine stop()
Definition: stop.f:20
subroutine preconvert2slapcol(irow, ia, jq, ja, nzs, nef)
Definition: preconvert2slapcol.f:26
subroutine dslugm(N, B, X, NELT, IA, JA, A, ISYM, NSAVE, ITOL, TOL, ITMAX, ITER, ERR, IERR, IUNIT, RWORK, LENW, IWORK, LENIW)
Definition: dslugm.f:4
subroutine extrapolate(yi, yn, ipkon, inum, kon, lakon, nfield, nk, ne, mi, ndim, orab, ielorien, co, iorienloc, cflag, vold, force, ielmat, thicke, ielprop, prop)
Definition: extrapolate.f:22
subroutine correctvel(hel, adv, vfa, ipnei, area, bv, xxn, neifa, lakonf, nef, neq)
Definition: correctvel.f:21
subroutine closefilefluid()
Definition: closefilefluid.f:20
subroutine calchcfa(nface, vfa, cocon, ncocon, ielmat, ntmat_, mi, ielfa, hcfa)
Definition: calchcfa.f:21
subroutine calcgamma(nface, ielfa, vel, gradvel, gamma, xlet, xxn, xxj, ipnei, betam, nef, flux)
Definition: calcgamma.f:21
subroutine calcstressheatflux(sti, umel, gradvel, qfx, hcel, gradtel, nef, isti, iqfx, mi)
Definition: calcstressheatflux.f:21
subroutine calcumfa(nface, vfa, shcon, nshcon, ielmat, ntmat_, ithermal, mi, ielfa, umfa)
Definition: calcumfa.f:21
subroutine predgmres(n, b, x, nelt, ia, ja, a, isym, itol, tol, itmax, iter, err, ierr, iunit, sb, sx, rgwk, lrgw, igwk, ligw, rwork, iwork)
Definition: predgmres.f:28
subroutine complete_hel_cyclic_blk(vel, hel, auv6, c, ipnei, neiel, neifa, ifatie, nef)
Definition: complete_hel_cyclic_blk.f:27
#define RENEW(a, b, c)
Definition: CalculiX.h:40
subroutine materialdata_cfd_comp(nef, vel, shcon, nshcon, ielmatf, ntmat_, mi, cvel, vfa, cocon, ncocon, physcon, cvfa, ithermal, nface, umel, umfa, ielfa, hcfa)
Definition: materialdata_cfd_comp.f:22
#define SFREE(a)
Definition: CalculiX.h:41
subroutine extrapol_pel(nface, ielfa, xrlfa, vel, vfa, ifabou, xbounact, nef, gradpel, gradpfa, neifa, rf, area, volume, xle, xxi, icyclic, xxn, ipnei, ifatie, coefmpc, nmpc, labmpc, ipompc, nodempc, ifaext, nfaext, nactdoh)
Definition: extrapol_pel.f:23
subroutine extrapol_kel(nface, ielfa, xrlfa, vel, vfa, ifabou, xbounact, nef, gradkel, gradkfa, neifa, rf, area, volume, xle, xxi, icyclic, xxn, ipnei, ifatie, xlet, xxj, coefmpc, nmpc, labmpc, ipompc, nodempc, ifaext, nfaext, nactdoh, umfa, physcon)
Definition: extrapol_kel.f:24
subroutine calcgammak(nface, ielfa, vel, gradkel, gamma, xlet, xxn, xxj, ipnei, betam, nef, flux)
Definition: calcgammak.f:21
subroutine calcrhoelcomp(nef, vel, shcon, ielmat, ntmat_, mi)
Definition: calcrhoelcomp.f:21
subroutine create_iau6(nef, ipnei, neiel, jq, irow, nzs, iau6, lakonf)
Definition: create_iau6.f:20
void mafillomain(ITG *nef, ITG *ipnei, ITG *neifa, ITG *neiel, double *vfa, double *xxn, double *area, double *au, double *ad, ITG *jq, ITG *irow, ITG *nzs, double *b, double *vel, double *umel, double *xlet, double *xle, double *gradofa, double *xxi, double *body, double *volume, ITG *ielfa, char *lakonf, ITG *ifabou, ITG *nbody, ITG *neq, double *dtimef, double *velo, double *veloo, double *cpfa, double *hcfa, double *cvel, double *gradvel, double *xload, double *gammat, double *xrlfa, double *xxj, ITG *nactdohinv, double *a1, double *a2, double *a3, double *flux, ITG *iau6, double *xxni, double *xxnj, ITG *iturbulent, double *gradkel, double *gradoel)
Definition: mafillomain.c:35
static double * f1
Definition: objectivemain_se.c:47
subroutine materialdata_cfd(nef, vel, shcon, nshcon, ielmatf, ntmat_, mi, cvel, vfa, cocon, ncocon, physcon, cvfa, ithermal, nface, umel, umfa, ielfa, hcfa, rhcon, nrhcon)
Definition: materialdata_cfd.f:22
subroutine extrapol_tel(nface, ielfa, xrlfa, vel, vfa, ifabou, xbounact, nef, gradtel, gradtfa, neifa, rf, area, volume, xle, xxi, icyclic, xxn, ipnei, ifatie, xload, xlet, xxj, coefmpc, nmpc, labmpc, ipompc, nodempc, ifaext, nfaext, nactdoh)
Definition: extrapol_tel.f:23
void mafillpmain(ITG *ne, char *lakonf, ITG *ipnei, ITG *neifa, ITG *neiel, double *vfa, double *area, double *adfa, double *xlet, double *cosa, double *volume, double *au, double *ad, ITG *jq, ITG *irow, double *ap, ITG *ielfa, ITG *ifabou, double *xle, double *b, double *xxn, ITG *neq, ITG *nzs, double *hfa, double *gradpel, double *bp, double *xxi, ITG *neij, double *xlen, double *cosb, ITG *iatleastonepressurebc, ITG *iau6, double *xxicn)
Definition: mafillpmain.c:33
void mafilltcompmain(ITG *nef, ITG *ipnei, ITG *neifa, ITG *neiel, double *vfa, double *xxn, double *area, double *au, double *ad, ITG *jq, ITG *irow, ITG *nzs, double *b, double *vel, double *umel, double *xlet, double *xle, double *gradtfa, double *xxi, double *body, double *volume, ITG *ielfa, char *lakonf, ITG *ifabou, ITG *nbody, ITG *neq, double *dtimef, double *velo, double *veloo, double *cpfa, double *hcfa, double *cvel, double *gradvel, double *xload, double *gammat, double *xrlfa, double *xxj, ITG *nactdohinv, double *a1, double *a2, double *a3, double *flux, ITG *iau6, double *xxni, double *xxnj)
Definition: mafilltcompmain.c:35
static ITG num_cpus
Definition: compfluid.c:37
subroutine complete_hel_cyclic(nef, bv, hel, adv, auv, jq, irow, ipnei, neiel, ifatie, c, lakonf, neifa, nzs)
Definition: complete_hel_cyclic.f:27
subroutine applyboun(ifaext, nfaext, ielfa, ikboun, ilboun, nboun, typeboun, nelemload, nload, sideload, isolidsurf, nsolidsurf, ifabou, nfabou, nface, nodeboun, ndirboun, ikmpc, ilmpc, labmpc, nmpc, nactdohinv, compressible, iatleastonepressurebc, ipkonf, kon, konf, nblk)
Definition: applyboun.f:24
subroutine printoutface(co, rhcon, nrhcon, ntmat_, vold, shcon, nshcon, cocon, ncocon, icompressible, istartset, iendset, ipkonf, lakonf, konf, ialset, prset, ttime, nset, set, nprint, prlab, ielmat, mi, ithermal, nactdoh, icfd, time, stn)
Definition: printoutface.f:23
void mastructf(ITG *nk, ITG *kon, ITG *ipkon, char *lakon, ITG *ne, ITG *icol, ITG *jq, ITG **mast1p, ITG **irowp, ITG *isolver, ITG *ipointer, ITG *nzs, ITG *ipnei, ITG *ineiel, ITG *mi)
Definition: mastructf.c:27
subroutine calcrhoel(nef, vel, rhcon, nrhcon, ielmat, ntmat_, ithermal, mi)
Definition: calcrhoel.f:21
subroutine correctvfa(nface, ielfa, area, vfa, ap, bp, xxn, ifabou, ipnei, nef, neifa, hfa, vel, xboun, lakonf, flux)
Definition: correctvfa.f:21
subroutine calcrhofacomp(nface, vfa, shcon, ielmat, ntmat_, mi, ielfa, ipnei, vel, nef, flux, gradpel, gradtel, xxj, betam, xlet)
Definition: calcrhofacomp.f:22
subroutine calccvfacomp(nface, vfa, shcon, nshcon, ielmat, ntmat_, mi, ielfa, cvfa, physcon)
Definition: calccvfacomp.f:21
void rhspmain(ITG *ne, char *lakon, ITG *ipnei, ITG *neifa, ITG *neiel, double *vfa, double *area, double *adfa, double *xlet, double *cosa, double *volume, double *au, double *ad, ITG *jq, ITG *irow, double *ap, ITG *ielfa, ITG *ifabou, double *xle, double *b, double *xxn, ITG *neq, ITG *nzs, double *hfa, double *gradpel, double *bp, double *xxi, ITG *neij, double *xlen, ITG *iatleastonepressurebc, double *xxicn)
Definition: rhspmain.c:34
void mafillvcompmain(ITG *nef, ITG *ipnei, ITG *neifa, ITG *neiel, double *vfa, double *xxn, double *area, double *au, double *ad, ITG *jq, ITG *irow, ITG *nzs, double *b, double *vel, double *cosa, double *umfa, double *xlet, double *xle, double *gradvfa, double *xxi, double *body, double *volume, ITG *ielfa, char *lakonf, ITG *ifabou, ITG *nbody, double *dtimef, double *velo, double *veloo, double *sel, double *xrlfa, double *gamma, double *xxj, ITG *nactdohinv, double *a1, double *a2, double *a3, double *flux, ITG *icyclic, double *c, ITG *ifatie, ITG *iau6, double *xxni, double *xxnj)
Definition: mafillvcompmain.c:35
subroutine calcgammao(nface, ielfa, vel, gradoel, gamma, xlet, xxn, xxj, ipnei, betam, nef, flux)
Definition: calcgammao.f:21
#define ITG
Definition: CalculiX.h:51
void mafillvmain(ITG *nef, ITG *ipnei, ITG *neifa, ITG *neiel, double *vfa, double *xxn, double *area, double *au, double *ad, ITG *jq, ITG *irow, ITG *nzs, double *b, double *vel, double *cosa, double *umfa, double *xlet, double *xle, double *gradvfa, double *xxi, double *body, double *volume, ITG *ielfa, char *lakonf, ITG *ifabou, ITG *nbody, double *dtimef, double *velo, double *veloo, double *sel, double *xrlfa, double *gamma, double *xxj, ITG *nactdohinv, double *a1, double *a2, double *a3, double *flux, ITG *icyclic, double *c, ITG *ifatie, ITG *iau6, double *xxni, double *xxnj, ITG *iturbulent, double *gradvel)
Definition: mafillvmain.c:35
subroutine initialcfd(nef, ipkonf, konf, lakonf, co, coel, cofa, nface, ielfa, area, ipnei, neiel, xxn, xxi, xle, xlen, xlet, xrlfa, cosa, volume, neifa, xxj, cosb, dmin, ifatie, cs, tieset, icyclic, c, neij, physcon, isolidsurf, nsolidsurf, dy, xxni, xxnj, xxicn, nflnei, iturbulent, rf)
Definition: initialcfd.f:24
#define NNEW(a, b, c)
Definition: CalculiX.h:39
subroutine fill_neiel(nef, ipnei, neiel, neielcp)
Definition: fill_neiel.f:20

Variable Documentation

◆ num_cpus

ITG num_cpus
static
Hosted by OpenAircraft.com, (Michigan UAV, LLC)