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

Go to the source code of this file.

Functions

void electromagnetics (double **cop, ITG *nk, ITG **konp, ITG **ipkonp, char **lakonp, ITG *ne, ITG *nodeboun, ITG *ndirboun, double *xboun, ITG *nboun, ITG **ipompcp, ITG **nodempcp, double **coefmpcp, char **labmpcp, ITG *nmpc, ITG *nodeforc, ITG *ndirforc, double *xforc, ITG *nforc, ITG **nelemloadp, char **sideloadp, double *xload, ITG *nload, ITG *nactdof, ITG **icolp, ITG *jq, ITG **irowp, ITG *neq, ITG *nzl, ITG *nmethod, ITG **ikmpcp, ITG **ilmpcp, ITG *ikboun, ITG *ilboun, double *elcon, ITG *nelcon, double *rhcon, ITG *nrhcon, double *alcon, ITG *nalcon, double *alzero, ITG **ielmatp, ITG **ielorienp, ITG *norien, double *orab, ITG *ntmat_, double *t0, double *t1, double *t1old, ITG *ithermal, double *prestr, ITG *iprestr, double **voldp, ITG *iperturb, double *sti, ITG *nzs, ITG *kode, char *filab, ITG *idrct, ITG *jmax, ITG *jout, double *timepar, double *eme, double *xbounold, double *xforcold, double *xloadold, double *veold, double *accold, char *amname, double *amta, ITG *namta, ITG *nam, ITG *iamforc, ITG **iamloadp, ITG *iamt1, double *alpha, ITG *iexpl, ITG *iamboun, double *plicon, ITG *nplicon, double *plkcon, ITG *nplkcon, double **xstatep, ITG *npmat_, ITG *istep, double *ttime, char *matname, double *qaold, ITG *mi, ITG *isolver, ITG *ncmat_, ITG *nstate_, ITG *iumat, double *cs, ITG *mcs, ITG *nkon, double **enerp, ITG *mpcinfo, char *output, double *shcon, ITG *nshcon, double *cocon, ITG *ncocon, double *physcon, ITG *nflow, double *ctrl, char **setp, ITG *nset, ITG **istartsetp, ITG **iendsetp, ITG **ialsetp, ITG *nprint, char *prlab, char *prset, ITG *nener, ITG *ikforc, ITG *ilforc, double *trab, ITG *inotr, ITG *ntrans, double **fmpcp, char *cbody, ITG *ibody, double *xbody, ITG *nbody, double *xbodyold, ITG *ielprop, double *prop, ITG *ntie, char **tiesetp, ITG *itpamp, ITG *iviewfile, char *jobnamec, double **tietolp, ITG *nslavs, double *thicke, ITG *ics, ITG *nalset, ITG *nmpc_, ITG *nmat, char *typeboun, ITG *iaxial, ITG *nload_, ITG *nprop, ITG *network, char *orname)
 

Function Documentation

◆ electromagnetics()

void electromagnetics ( double **  cop,
ITG nk,
ITG **  konp,
ITG **  ipkonp,
char **  lakonp,
ITG ne,
ITG nodeboun,
ITG ndirboun,
double *  xboun,
ITG nboun,
ITG **  ipompcp,
ITG **  nodempcp,
double **  coefmpcp,
char **  labmpcp,
ITG nmpc,
ITG nodeforc,
ITG ndirforc,
double *  xforc,
ITG nforc,
ITG **  nelemloadp,
char **  sideloadp,
double *  xload,
ITG nload,
ITG nactdof,
ITG **  icolp,
ITG jq,
ITG **  irowp,
ITG neq,
ITG nzl,
ITG nmethod,
ITG **  ikmpcp,
ITG **  ilmpcp,
ITG ikboun,
ITG ilboun,
double *  elcon,
ITG nelcon,
double *  rhcon,
ITG nrhcon,
double *  alcon,
ITG nalcon,
double *  alzero,
ITG **  ielmatp,
ITG **  ielorienp,
ITG norien,
double *  orab,
ITG ntmat_,
double *  t0,
double *  t1,
double *  t1old,
ITG ithermal,
double *  prestr,
ITG iprestr,
double **  voldp,
ITG iperturb,
double *  sti,
ITG nzs,
ITG kode,
char *  filab,
ITG idrct,
ITG jmax,
ITG jout,
double *  timepar,
double *  eme,
double *  xbounold,
double *  xforcold,
double *  xloadold,
double *  veold,
double *  accold,
char *  amname,
double *  amta,
ITG namta,
ITG nam,
ITG iamforc,
ITG **  iamloadp,
ITG iamt1,
double *  alpha,
ITG iexpl,
ITG iamboun,
double *  plicon,
ITG nplicon,
double *  plkcon,
ITG nplkcon,
double **  xstatep,
ITG npmat_,
ITG istep,
double *  ttime,
char *  matname,
double *  qaold,
ITG mi,
ITG isolver,
ITG ncmat_,
ITG nstate_,
ITG iumat,
double *  cs,
ITG mcs,
ITG nkon,
double **  enerp,
ITG mpcinfo,
char *  output,
double *  shcon,
ITG nshcon,
double *  cocon,
ITG ncocon,
double *  physcon,
ITG nflow,
double *  ctrl,
char **  setp,
ITG nset,
ITG **  istartsetp,
ITG **  iendsetp,
ITG **  ialsetp,
ITG nprint,
char *  prlab,
char *  prset,
ITG nener,
ITG ikforc,
ITG ilforc,
double *  trab,
ITG inotr,
ITG ntrans,
double **  fmpcp,
char *  cbody,
ITG ibody,
double *  xbody,
ITG nbody,
double *  xbodyold,
ITG ielprop,
double *  prop,
ITG ntie,
char **  tiesetp,
ITG itpamp,
ITG iviewfile,
char *  jobnamec,
double **  tietolp,
ITG nslavs,
double *  thicke,
ITG ics,
ITG nalset,
ITG nmpc_,
ITG nmat,
char *  typeboun,
ITG iaxial,
ITG nload_,
ITG nprop,
ITG network,
char *  orname 
)
79  {
80 
81  char description[13]=" ",*lakon=NULL,jobnamef[396]="",
82  *labmpc=NULL,kind1[2]="E",kind2[2]="E",*set=NULL,*tieset=NULL,
83  cflag[1]=" ",*sideloadref=NULL,*sideload=NULL;
84 
85  ITG *inum=NULL,k,iout=0,icntrl,iinc=0,jprint=0,iit=-1,jnz=0,
86  icutb=0,istab=0,ifreebody,uncoupled=0,maxfaces,indexe,nope,
87  iperturb_sav[2],*icol=NULL,*irow=NULL,ielas=0,icmd=0,j,
88  memmpc_,mpcfree,icascade,maxlenmpc,*nodempc=NULL,*iaux=NULL,
89  *itg=NULL,*ineighe=NULL,null=0,iactive[3],neqterms,ntflag,
90  *ieg=NULL,ntg=0,ntr,*kontri=NULL,*nloadtr=NULL,index,
91  *ipiv=NULL,ntri,newstep,mode=-1,noddiam=-1,nasym=0,
92  ntrit,*inocs=NULL,inewton=0,*ipobody=NULL,*nacteq=NULL,
93  *nactdog=NULL,nteq,nmastnode,imast,massact[2],
94  *ipkon=NULL,*kon=NULL,*ielorien=NULL,nmethodact,ne2=0,
95  *ielmat=NULL,inext,itp=0,symmetryflag=0,inputformat=0,
96  iitterm=0,ngraph=1,ithermalact=2,*islavact=NULL,neini,
97  *ipompc=NULL,*ikmpc=NULL,*ilmpc=NULL,i0ref,irref,icref,
98  *islavnode=NULL,*imastnode=NULL,*nslavnode=NULL,mortar=0,
99  mt=mi[1]+1,*nactdofinv=NULL, inode,idir,*islavsurf=NULL,
100  iemchange=0,nzsrad,*mast1rad=NULL,*irowrad=NULL,*icolrad=NULL,
101  *jqrad=NULL,*ipointerrad=NULL,*integerglob=NULL,im,ne0,
102  mass[2]={0,0}, stiffness=1, buckling=0, rhsi=1, intscheme=0,idiscon=0,
103  coriolis=0,*ipneigh=NULL,*neigh=NULL,i,icfd=0,id,node,networknode,
104  iflagact=0,*nodorig=NULL,*ipivr=NULL,*inomat=NULL,*nodface=NULL,
105  *ipoface=NULL,*istartset=NULL,*iendset=NULL,*ialset=NULL,
106  *nelemloadref=NULL,*iamloadref=NULL,nloadref,kscale=1,
107  *nelemload=NULL,*iamload=NULL,*idefload=NULL,ialeatoric=0,
108  *iponoel=NULL,*inoel=NULL,inoelsize;
109 
110  double *stn=NULL,*v=NULL,*een=NULL,cam[5],*epn=NULL,*cdn=NULL,
111  *f=NULL,*fn=NULL,qa[4]={0.,0.,-1.,0.},qam[2]={0.,0.},dtheta,theta,
112  err,ram[4]={0.,0.,0.,0.},*springarea=NULL,*h0=NULL,
113  ram1[2]={0.,0.},ram2[2]={0.,0.},deltmx,*clearini=NULL,
114  uam[2]={0.,0.},*vini=NULL,*ac=NULL,qa0,qau,ea,ptime,
115  *t1act=NULL,qamold[2],*xbounact=NULL,*bc=NULL,
116  *xforcact=NULL,*xloadact=NULL,*fext=NULL,h0scale=1.,
117  reltime,time,bet=0.,gam=0.,*aux1=NULL,*aux2=NULL,dtime,*fini=NULL,
118  *fextini=NULL,*veini=NULL,*xstateini=NULL,*h0ref=NULL,
119  *ampli=NULL,*eei=NULL,*t1ini=NULL,*tinc,*tper,*tmin,*tmax,
120  *xbounini=NULL,*xstiff=NULL,*stx=NULL,*cv=NULL,*cvini=NULL,
121  *enern=NULL,*coefmpc=NULL,*aux=NULL,*xstaten=NULL,
122  *enerini=NULL,*emn=NULL,*xmastnor=NULL,*fnext=NULL,
123  *tarea=NULL,*tenv=NULL,*erad=NULL,*fnr=NULL,*fni=NULL,
124  *adview=NULL,*auview=NULL,*qfx=NULL,*adaux=NULL,
125  *qfn=NULL,*co=NULL,*vold=NULL,*fenv=NULL,sigma=0.,
126  *xbodyact=NULL,*cgr=NULL,dthetaref,*vr=NULL,*vi=NULL,
127  *stnr=NULL,*stni=NULL,*vmax=NULL,*stnmax=NULL,*fmpc=NULL,*ener=NULL,
128  *f_cm=NULL,*f_cs=NULL,*tietol=NULL,
129  *xstate=NULL,*eenmax=NULL,*adrad=NULL,*aurad=NULL,*bcr=NULL,
130  *emeini=NULL,*doubleglob=NULL,*au=NULL,
131  *ad=NULL,*b=NULL,*aub=NULL,*adb=NULL,*pslavsurf=NULL,*pmastsurf=NULL,
132  *cdnr=NULL,*cdni=NULL,*energyini=NULL,*energy=NULL;
133 
134 #ifdef SGI
135  ITG token;
136 #endif
137 
138  ne0=*ne;
139 
140  /* next line is needed to avoid that elements with negative ipkon
141  are taken into account in extrapolate.f */
142 
143  strcpy1(&filab[2],"C",1);
144 
145  if(*nmethod==8){
146  *nmethod=1;
147  }else if(*nmethod==9){
148  *nmethod=4;
149  }else if(*nmethod==10){
150  *nmethod=2;
151  }
152 
153  for(k=0;k<3;k++){
154  strcpy1(&jobnamef[k*132],&jobnamec[k*132],132);
155  }
156 
157  qa0=ctrl[20];qau=ctrl[21];ea=ctrl[23];deltmx=ctrl[26];
158  i0ref=ctrl[0];irref=ctrl[1];icref=ctrl[3];
159 
160  memmpc_=mpcinfo[0];mpcfree=mpcinfo[1];icascade=mpcinfo[2];
161  maxlenmpc=mpcinfo[3];
162 
163  icol=*icolp;irow=*irowp;co=*cop;vold=*voldp;
164  ipkon=*ipkonp;lakon=*lakonp;kon=*konp;ielorien=*ielorienp;
165  ielmat=*ielmatp;ener=*enerp;xstate=*xstatep;
166 
167  ipompc=*ipompcp;labmpc=*labmpcp;ikmpc=*ikmpcp;ilmpc=*ilmpcp;
168  fmpc=*fmpcp;nodempc=*nodempcp;coefmpc=*coefmpcp;
169 
170  set=*setp;istartset=*istartsetp;iendset=*iendsetp;ialset=*ialsetp;
171  tieset=*tiesetp;tietol=*tietolp;
172 
173  nelemload=*nelemloadp;iamload=*iamloadp;
174  sideload=*sideloadp;
175 
176  tinc=&timepar[0];
177  tper=&timepar[1];
178  tmin=&timepar[2];
179  tmax=&timepar[3];
180 
181  /* invert nactdof */
182 
183  NNEW(nactdofinv,ITG,mt**nk);
184  NNEW(nodorig,ITG,*nk);
185  FORTRAN(gennactdofinv,(nactdof,nactdofinv,nk,mi,nodorig,
186  ipkon,lakon,kon,ne));
187  SFREE(nodorig);
188 
189  /* allocating a field for the stiffness matrix */
190 
191  NNEW(xstiff,double,(long long)27*mi[0]**ne);
192 
193  /* allocating force fields */
194 
195  NNEW(f,double,neq[1]);
196  NNEW(fext,double,neq[1]);
197 
198  NNEW(b,double,neq[1]);
199  NNEW(vini,double,mt**nk);
200 
201  NNEW(aux,double,7*maxlenmpc);
202  NNEW(iaux,ITG,maxlenmpc);
203 
204  /* allocating fields for the actual external loading */
205 
206  NNEW(xbounact,double,*nboun);
207  NNEW(xbounini,double,*nboun);
208  for(k=0;k<*nboun;++k){xbounact[k]=xbounold[k];}
209  NNEW(xforcact,double,*nforc);
210  NNEW(xloadact,double,2**nload);
211  NNEW(xbodyact,double,7**nbody);
212  /* copying the rotation axis and/or acceleration vector */
213  for(k=0;k<7**nbody;k++){xbodyact[k]=xbody[k];}
214 
215  /* assigning the body forces to the elements */
216 
217  if(*nbody>0){
218  ifreebody=*ne+1;
219  NNEW(ipobody,ITG,2*ifreebody**nbody);
220  for(k=1;k<=*nbody;k++){
221  FORTRAN(bodyforce,(cbody,ibody,ipobody,nbody,set,istartset,
222  iendset,ialset,&inewton,nset,&ifreebody,&k));
223  RENEW(ipobody,ITG,2*(*ne+ifreebody));
224  }
225  RENEW(ipobody,ITG,2*(ifreebody-1));
226  }
227 
228  /* for thermal calculations: forced convection and cavity
229  radiation*/
230 
231  if(*ithermal>1){
232  NNEW(itg,ITG,*nload+3**nflow);
233  NNEW(ieg,ITG,*nflow);
234  /* max 6 triangles per face, 4 entries per triangle */
235  NNEW(kontri,ITG,24**nload);
236  NNEW(nloadtr,ITG,*nload);
237  NNEW(nacteq,ITG,4**nk);
238  NNEW(nactdog,ITG,4**nk);
239  NNEW(v,double,mt**nk);
240  FORTRAN(envtemp,(itg,ieg,&ntg,&ntr,sideload,nelemload,
241  ipkon,kon,lakon,ielmat,ne,nload,
242  kontri,&ntri,nloadtr,nflow,ndirboun,nactdog,
243  nodeboun,nacteq,nboun,ielprop,prop,&nteq,
244  v,network,physcon,shcon,ntmat_,co,
245  vold,set,nshcon,rhcon,nrhcon,mi,nmpc,nodempc,
246  ipompc,labmpc,ikboun,&nasym,ttime,&time,iaxial));
247  SFREE(v);
248 
249  if((*mcs>0)&&(ntr>0)){
250  NNEW(inocs,ITG,*nk);
251  radcyc(nk,kon,ipkon,lakon,ne,cs,mcs,nkon,ialset,istartset,
252  iendset,&kontri,&ntri,&co,&vold,&ntrit,inocs,mi);
253  }
254  else{ntrit=ntri;}
255 
256  nzsrad=100*ntr;
257  NNEW(mast1rad,ITG,nzsrad);
258  NNEW(irowrad,ITG,nzsrad);
259  NNEW(icolrad,ITG,ntr);
260  NNEW(jqrad,ITG,ntr+1);
261  NNEW(ipointerrad,ITG,ntr);
262 
263  if(ntr>0){
264  mastructrad(&ntr,nloadtr,sideload,ipointerrad,
265  &mast1rad,&irowrad,&nzsrad,
266  jqrad,icolrad);
267  }
268 
269  /* determine the network elements belonging to a given node (for usage
270  in user subroutine film */
271 
272 // if(ntg>0){
273  if((*network>0)||(ntg>0)){
274  NNEW(iponoel,ITG,*nk);
275  NNEW(inoel,ITG,2**nkon);
276  if(*network>0){
277  FORTRAN(networkelementpernode,(iponoel,inoel,lakon,ipkon,kon,
278  &inoelsize,nflow,ieg,ne,network));
279  }
280  RENEW(inoel,ITG,2*inoelsize);
281  }
282 
283  SFREE(ipointerrad);SFREE(mast1rad);
284  RENEW(irowrad,ITG,nzsrad);
285 
286  RENEW(itg,ITG,ntg);
287  NNEW(ineighe,ITG,ntg);
288  RENEW(kontri,ITG,4*ntrit);
289  RENEW(nloadtr,ITG,ntr);
290 
291  NNEW(adview,double,ntr);
292  NNEW(auview,double,2*nzsrad);
293  NNEW(tarea,double,ntr);
294  NNEW(tenv,double,ntr);
295  NNEW(fenv,double,ntr);
296  NNEW(erad,double,ntr);
297 
298  NNEW(ac,double,nteq*nteq);
299  NNEW(bc,double,nteq);
300  NNEW(ipiv,ITG,nteq);
301  NNEW(adrad,double,ntr);
302  NNEW(aurad,double,2*nzsrad);
303  NNEW(bcr,double,ntr);
304  NNEW(ipivr,ITG,ntr);
305  }
306  if(*ithermal>1){NNEW(qfx,double,3*mi[0]**ne);}
307 
308  /* allocating a field for the instantaneous amplitude */
309 
310  NNEW(ampli,double,*nam);
311 
312  NNEW(fini,double,neq[1]);
313 
314  /* allocating fields for nonlinear dynamics */
315 
316  if(*nmethod==4){
317  mass[0]=1;
318  mass[1]=1;
319  NNEW(aux2,double,neq[1]);
320  NNEW(fextini,double,neq[1]);
321  NNEW(cv,double,neq[1]);
322  NNEW(cvini,double,neq[1]);
323  NNEW(veini,double,mt**nk);
324  NNEW(adb,double,neq[1]);
325  NNEW(aub,double,nzs[1]);
326  }
327 
328  qa[0]=qaold[0];
329  qa[1]=qaold[1];
330 
331  /* normalizing the time */
332 
333  FORTRAN(checktime,(itpamp,namta,tinc,ttime,amta,tmin,&inext,&itp,istep,tper));
334  dtheta=(*tinc)/(*tper);
335  dthetaref=dtheta;
336  if(dtheta<=1.e-6){
337  printf("\n *ERROR in nonlingeo\n");
338  printf(" increment size smaller than one millionth of step size\n");
339  printf(" increase increment size\n\n");
340  }
341  *tmin=*tmin/(*tper);
342  *tmax=*tmax/(*tper);
343  theta=0.;
344 
345  /* calculating an initial flux norm */
346 
347  if(*ithermal!=2){
348  if(qau>1.e-10){qam[0]=qau;}
349  else if(qa0>1.e-10){qam[0]=qa0;}
350  else if(qa[0]>1.e-10){qam[0]=qa[0];}
351  else {qam[0]=1.e-2;}
352  }
353  if(*ithermal>1){
354  if(qau>1.e-10){qam[1]=qau;}
355  else if(qa0>1.e-10){qam[1]=qa0;}
356  else if(qa[1]>1.e-10){qam[1]=qa[1];}
357  else {qam[1]=1.e-2;}
358  }
359 
360 
361  /*********************************************************************/
362 
363  /* calculating the initial quasi-static magnetic intensity due to
364  the coil current */
365 
366  /*********************************************************************/
367 
368  /* calculate the current density in the coils
369 
370  in this section nload, nforc, nbody and nam are set to zero; the
371  electrical potential is supposed to be given (in the form of a
372  temperature), the current is calculated (in the form of heat
373  flux) by thermal analogy */
374 
375  reltime=1.;
376  time=0.;
377  dtime=0.;
378  ithermalact=2;
379 
380  nmethodact=1;
381  massact[0]=0;
382  massact[1]=0;
383 
384  if(*ithermal<=1){
385  NNEW(qfx,double,3*mi[0]**ne);
386  NNEW(t0,double,*nk);
387  }
388  if(strcmp1(&filab[3567],"ECD ")==0){NNEW(qfn,double,3**nk);}
389 
390  /* the coil current is assumed to be applied at once, i.e. as
391  step loading; the calculation, however, is a quasi-static
392  calculation */
393 
394  FORTRAN(tempload,(xforcold,xforc,xforcact,iamforc,&null,xloadold,xload,
395  xloadact,iamload,&null,ibody,xbody,&null,xbodyold,xbodyact,
396  t1old,t1,t1act,iamt1,nk,amta,
397  namta,&null,ampli,&time,&reltime,ttime,&dtime,&ithermalact,nmethod,
398  xbounold,xboun,xbounact,iamboun,nboun,
399  nodeboun,ndirboun,nodeforc,ndirforc,istep,&iinc,
400  co,vold,itg,&ntg,amname,ikboun,ilboun,nelemload,sideload,mi,
401  ntrans,trab,inotr,veold,integerglob,doubleglob,tieset,istartset,
402  iendset,ialset,ntie,nmpc,ipompc,ikmpc,ilmpc,nodempc,coefmpc,
403  ipobody,iponoel,inoel));
404 
405  cam[0]=0.;cam[1]=0.;
406 
407  /* deactivating all elements except the shells */
408 
409  for(i=0;i<*ne;i++){
410  if(strcmp1(&lakon[8*i+6],"L")!=0){
411  ipkon[i]=-ipkon[i]-2;
412  }
413  }
414 
415  remastruct(ipompc,&coefmpc,&nodempc,nmpc,
416  &mpcfree,nodeboun,ndirboun,nboun,ikmpc,ilmpc,ikboun,ilboun,
417  labmpc,nk,&memmpc_,&icascade,&maxlenmpc,
418  kon,ipkon,lakon,ne,nactdof,icol,jq,&irow,isolver,
419  neq,nzs,&nmethodact,&f,&fext,&b,&aux2,&fini,&fextini,
420  &adb,&aub,&ithermalact,iperturb,mass,mi,iexpl,&mortar,
421  typeboun,&cv,&cvini,&iit,network);
422 
423  /* invert nactdof */
424 
425  SFREE(nactdofinv);
426  NNEW(nactdofinv,ITG,mt**nk);
427  NNEW(nodorig,ITG,*nk);
428  FORTRAN(gennactdofinv,(nactdof,nactdofinv,nk,mi,nodorig,
429  ipkon,lakon,kon,ne));
430  SFREE(nodorig);
431 
432  iout=-1;
433 
434  NNEW(fn,double,mt**nk);
435  NNEW(inum,ITG,*nk);
436  NNEW(v,double,mt**nk);
437 
438  results(co,nk,kon,ipkon,lakon,ne,v,stn,inum,stx,
439  elcon,nelcon,rhcon,nrhcon,alcon,nalcon,alzero,ielmat,
440  ielorien,norien,orab,ntmat_,t0,t1act,&ithermalact,
441  prestr,iprestr,filab,eme,emn,een,iperturb,f,fn,nactdof,&iout,
442  qa,vold,b,nodeboun,ndirboun,xbounact,nboun,ipompc,nodempc,coefmpc,
443  labmpc,nmpc,&nmethodact,cam,&neq[1],veold,accold,&bet,
444  &gam,&dtime,&time,ttime,plicon,nplicon,plkcon,nplkcon,
445  xstateini,xstiff,xstate,npmat_,epn,matname,mi,&ielas,&icmd,
446  ncmat_,nstate_,sti,vini,ikboun,ilboun,ener,enern,emeini,xstaten,
447  eei,enerini,alcon,nalcon,set,nset,istartset,iendset,
448  ialset,nprint,prlab,prset,qfx,qfn,trab,inotr,ntrans,fmpc,
449  nelemload,&null,ikmpc,ilmpc,istep,&iinc,springarea,&reltime,
450  ne,xforc,&null,thicke,shcon,nshcon,
451  sideload,xloadact,xloadold,&icfd,inomat,pslavsurf,pmastsurf,
452  &mortar,islavact,cdn,islavnode,nslavnode,ntie,clearini,
453  islavsurf,ielprop,prop,energyini,energy,&kscale,iponoel,
454  inoel,nener,orname,network,ipobody,xbodyact,ibody);
455 
456  SFREE(fn);SFREE(inum);SFREE(v);
457 
458  iout=1;
459 
460  NNEW(ad,double,neq[1]);
461  NNEW(au,double,nzs[1]);
462 
463  mafillsmmain(co,nk,kon,ipkon,lakon,ne,nodeboun,ndirboun,xbounold,nboun,
464  ipompc,nodempc,coefmpc,nmpc,nodeforc,ndirforc,xforcact,
465  &null,nelemload,sideload,xloadact,&null,xbodyact,ipobody,
466  &null,cgr,ad,au,fext,nactdof,icol,jq,irow,neq,nzl,
467  &nmethodact,ikmpc,ilmpc,ikboun,ilboun,
468  elcon,nelcon,rhcon,nrhcon,alcon,nalcon,alzero,
469  ielmat,ielorien,norien,orab,ntmat_,
470  t0,t1act,&ithermalact,prestr,iprestr,vold,iperturb,sti,
471  nzs,stx,adb,aub,iexpl,plicon,nplicon,plkcon,nplkcon,
472  xstiff,npmat_,&dtime,matname,mi,
473  ncmat_,massact,&stiffness,&buckling,&rhsi,&intscheme,
474  physcon,shcon,nshcon,alcon,nalcon,ttime,&time,istep,&iinc,
475  &coriolis,ibody,xloadold,&reltime,veold,springarea,nstate_,
476  xstateini,xstate,thicke,integerglob,doubleglob,
477  tieset,istartset,iendset,ialset,ntie,&nasym,pslavsurf,
478  pmastsurf,&mortar,clearini,ielprop,prop,&ne0,fnext,&kscale,
479  iponoel,inoel,network);
480 
481  if(nmethodact==0){
482 
483  /* error occurred in mafill: storing the geometry in frd format */
484 
485  ++*kode;
486 
487  ptime=*ttime+time;
488  frd(co,nk,kon,ipkon,lakon,ne,v,stn,inum,&nmethodact,
489  kode,filab,een,t1,fn,&ptime,epn,ielmat,matname,enern,xstaten,
490  nstate_,istep,&iinc,&ithermalact,qfn,&mode,&noddiam,trab,inotr,
491  ntrans,orab,ielorien,norien,description,ipneigh,neigh,
492  mi,sti,vr,vi,stnr,stni,vmax,stnmax,&ngraph,veold,ener,ne,
493  cs,set,nset,istartset,iendset,ialset,eenmax,fnr,fni,emn,
494  thicke,jobnamec,output,qfx,cdn,&mortar,cdnr,cdni,nmat);
495 
496  FORTRAN(stop,());
497 
498  }
499 
500  for(k=0;k<neq[1];++k){
501  b[k]=fext[k]-f[k];
502  }
503 
504  if(*isolver==0){
505 #ifdef SPOOLES
506  spooles(ad,au,adb,aub,&sigma,b,icol,irow,&neq[1],&nzs[1],
507  &symmetryflag,&inputformat,&nzs[2]);
508 #else
509  printf("*ERROR in nonlingeo: the SPOOLES library is not linked\n\n");
510  FORTRAN(stop,());
511 #endif
512  }
513  else if((*isolver==2)||(*isolver==3)){
514  preiter(ad,&au,b,&icol,&irow,&neq[1],&nzs[1],isolver,iperturb);
515  }
516  else if(*isolver==4){
517 #ifdef SGI
518  token=1;
519  sgi_main(ad,au,adb,aub,&sigma,b,icol,irow,&neq[1],&nzs[1],token);
520 #else
521  printf("*ERROR in nonlingeo: the SGI library is not linked\n\n");
522  FORTRAN(stop,());
523 #endif
524  }
525  else if(*isolver==5){
526 #ifdef TAUCS
527  tau(ad,&au,adb,aub,&sigma,b,icol,&irow,&neq[1],&nzs[1]);
528 #else
529  printf("*ERROR in nonlingeo: the TAUCS library is not linked\n\n");
530  FORTRAN(stop,());
531 #endif
532  }
533  else if(*isolver==7){
534 #ifdef PARDISO
535  pardiso_main(ad,au,adb,aub,&sigma,b,icol,irow,&neq[1],&nzs[1],
536  &symmetryflag,&inputformat,jq,&nzs[2]);
537 #else
538  printf("*ERROR in nonlingeo: the PARDISO library is not linked\n\n");
539  FORTRAN(stop,());
540 #endif
541  }
542 
543  SFREE(au);SFREE(ad);
544 
545  NNEW(v,double,mt**nk);
546 // memcpy(&v[0],&vold[0],sizeof(double)*mt**nk);
547 
548  NNEW(fn,double,mt**nk);
549 
550  NNEW(inum,ITG,*nk);
551  results(co,nk,kon,ipkon,lakon,ne,v,stn,inum,stx,
552  elcon,nelcon,rhcon,nrhcon,alcon,nalcon,alzero,ielmat,
553  ielorien,norien,orab,ntmat_,t0,t1act,&ithermalact,
554  prestr,iprestr,filab,eme,emn,een,iperturb,
555  f,fn,nactdof,&iout,qa,vold,b,nodeboun,
556  ndirboun,xbounact,nboun,ipompc,
557  nodempc,coefmpc,labmpc,nmpc,&nmethodact,cam,&neq[1],veold,accold,
558  &bet,&gam,&dtime,&time,ttime,plicon,nplicon,plkcon,nplkcon,
559  xstateini,xstiff,xstate,npmat_,epn,matname,mi,&ielas,
560  &icmd,ncmat_,nstate_,sti,vini,ikboun,ilboun,ener,enern,
561  emeini,xstaten,eei,enerini,alcon,nalcon,set,nset,istartset,
562  iendset,ialset,nprint,prlab,prset,qfx,qfn,trab,inotr,ntrans,
563  fmpc,nelemload,&null,ikmpc,ilmpc,istep,&iinc,springarea,
564  &reltime,ne,xforc,&null,thicke,shcon,nshcon,
565  sideload,xloadact,xloadold,&icfd,inomat,pslavsurf,pmastsurf,
566  &mortar,islavact,cdn,islavnode,nslavnode,ntie,clearini,
567  islavsurf,ielprop,prop,energyini,energy,&kscale,iponoel,
568  inoel,nener,orname,network,ipobody,xbodyact,ibody);
569 
570 // memcpy(&vold[0],&v[0],sizeof(double)*mt**nk);
571 
572  /* reactivating the non-shell elements (for mesh output purposes)
573  deactivating the initial temperature for the non-shell nodes */
574 
575  for(i=0;i<*ne;i++){
576  if(strcmp1(&lakon[8*i+6],"L")!=0){
577  ipkon[i]=-ipkon[i]-2;
578  }else if(ipkon[i]!=-1){
579 
580  /* copy shell results */
581 
582  indexe=ipkon[i];
583  if(strcmp1(&lakon[8*i+3],"6")==0){nope=6;}
584  else if(strcmp1(&lakon[8*i+3],"8")==0){nope=8;}
585  else if(strcmp1(&lakon[8*i+3],"1")==0){nope=15;}
586  else{nope=20;}
587  for(j=0;j<nope;j++){
588  node=kon[indexe+j];
589  vold[mt*(node-1)]=v[mt*(node-1)];
590  }
591  }
592  }
593 
594  /* deactivating the output of temperatures */
595 
596  if(strcmp1(&filab[87],"NT ")==0){
597  ntflag=1;
598  strcpy1(&filab[87]," ",4);
599  }else{ntflag=0;}
600 
601  /* createinum is called in order to store the nodes and elements
602  of the complete structure, not only of the coil */
603 
604  FORTRAN(createinum,(ipkon,inum,kon,lakon,nk,ne,&cflag[0],nelemload,
605  nload,nodeboun,nboun,ndirboun,ithermal,co,vold,mi,ielmat));
606 
607  ++*kode;
608  if(*mcs!=0){
609  ptime=*ttime+time;
610  frdcyc(co,nk,kon,ipkon,lakon,ne,v,stn,inum,&nmethodact,kode,filab,een,
611  t1act,fn,&ptime,epn,ielmat,matname,cs,mcs,nkon,enern,xstaten,
612  nstate_,istep,&iinc,iperturb,ener,mi,output,&ithermalact,qfn,
613  ialset,istartset,iendset,trab,inotr,ntrans,orab,ielorien,
614  norien,stx,veold,&noddiam,set,nset,emn,thicke,jobnamec,ne,
615  cdn,&mortar,nmat,qfx);
616  }else{
617 
618  ptime=*ttime+time;
619  frd(co,nk,kon,ipkon,lakon,ne,v,stn,inum,&nmethodact,
620  kode,filab,een,t1act,fn,&ptime,epn,ielmat,matname,enern,xstaten,
621  nstate_,istep,&iinc,&ithermalact,qfn,&mode,&noddiam,trab,inotr,
622  ntrans,orab,ielorien,norien,description,ipneigh,neigh,
623  mi,stx,vr,vi,stnr,stni,vmax,stnmax,&ngraph,veold,ener,ne,
624  cs,set,nset,istartset,iendset,ialset,eenmax,fnr,fni,emn,
625  thicke,jobnamec,output,qfx,cdn,&mortar,cdnr,cdni,nmat);
626 
627  }
628  SFREE(inum);SFREE(v);SFREE(fn);
629 
630  /* reactivating the temperature output, if previously deactivated */
631 
632  if(ntflag==1){
633  strcpy1(&filab[87],"NT ",4);
634  }
635 
636  NNEW(inomat,ITG,*nk);
637 
638  /* calculating the magnetic intensity caused by the current */
639 
640  FORTRAN(assigndomtonodes,(ne,lakon,ipkon,kon,ielmat,inomat,
641  elcon,ncmat_,ntmat_,mi,&ne2));
642 
643  NNEW(h0ref,double,3**nk);
644  NNEW(h0,double,3**nk);
645 
646  biosav(ipkon,kon,lakon,ne,co,qfx,h0ref,mi,inomat,nk);
647 
648  if(*ithermal<=1){SFREE(qfx);SFREE(t0);}
649  if(strcmp1(&filab[3567],"ECD ")==0)SFREE(qfn);
650 
651  /* deactivating the shell elements */
652 
653  for(i=0;i<*ne;i++){
654  if(strcmp1(&lakon[8*i+6],"L")==0){
655  ipkon[i]=-ipkon[i]-2;
656  }
657  }
658 
659 /**************************************************************/
660 /* creating connecting MPC's between the domains */
661 /**************************************************************/
662 
663 /* creating contact ties between the domains */
664 
665  if(*istep==1){
666 
667  NNEW(nodface,ITG,5*6**ne);
668  NNEW(ipoface,ITG,*nk);
669 
670  RENEW(set,char,81*(*nset+3));
671  RENEW(istartset,ITG,*nset+3);
672  RENEW(iendset,ITG,*nset+3);
673  RENEW(ialset,ITG,*nalset+6**ne);
674  RENEW(tieset,char,243*(*ntie+5));
675  RENEW(tietol,double,3*(*ntie+5));
676 
677  FORTRAN(createtiedsurfs,(nodface,ipoface,set,istartset,
678  iendset,ialset,tieset,inomat,ne,ipkon,lakon,kon,ntie,
679  tietol,nalset,nk,nset,iactive));
680 
681  SFREE(nodface);SFREE(ipoface);
682  RENEW(set,char,81**nset);
683  RENEW(istartset,ITG,*nset);
684  RENEW(iendset,ITG,*nset);
685  RENEW(ialset,ITG,*nalset);
686  RENEW(tieset,char,243**ntie);
687  RENEW(tietol,double,3**ntie);
688 
689  /* tied contact constraints: generate appropriate MPC's */
690 
691  tiedcontact(ntie,tieset,nset,set,istartset,iendset,ialset,
692  lakon,ipkon,kon,tietol,nmpc, &mpcfree,&memmpc_,
693  &ipompc,&labmpc,&ikmpc,&ilmpc,&fmpc,&nodempc,&coefmpc,
694  ithermal,co,vold,&icfd,nmpc_,mi,nk,istep,ikboun,nboun,
695  kind1,kind2);
696 
697  /* mapping h0ref from the phi domain onto the border of
698  the A and A-V domains */
699 
700  FORTRAN(calch0interface,(nmpc,ipompc,nodempc,coefmpc,h0ref));
701 
702  }
703 
704 /**************************************************************/
705 /* creating the A.n MPC */
706 /**************************************************************/
707 
708  /* identifying the interfaces between the A and A-V domains
709  and the phi-domain */
710 
711  FORTRAN(generateeminterfaces,(istartset,iendset,
712  ialset,iactive,ipkon,lakon,kon,ikmpc,nmpc,&maxfaces));
713 
714  for(i=1;i<3;i++){
715  imast=iactive[i];
716  if(imast==0) continue;
717 
718  /* determining the normals on the face */
719 
720  NNEW(imastnode,ITG,8*maxfaces);
721  NNEW(xmastnor,double,3*8*maxfaces);
722 
723  FORTRAN(normalsoninterface,(istartset,iendset,
724  ialset,&imast,ipkon,kon,lakon,imastnode,&nmastnode,
725  xmastnor,co));
726 
727  /* enlarging the fields for MPC's */
728 
729  *nmpc_=*nmpc_+nmastnode;
730  RENEW(ipompc,ITG,*nmpc_);
731  RENEW(labmpc,char,20**nmpc_+1);
732  RENEW(ikmpc,ITG,*nmpc_);
733  RENEW(ilmpc,ITG,*nmpc_);
734  RENEW(fmpc,double,*nmpc_);
735 
736  /* determining the maximum number of terms;
737  expanding nodempc and coefmpc to accommodate
738  those terms */
739 
740  neqterms=3*nmastnode;
741  index=memmpc_;
742  (memmpc_)+=neqterms;
743  RENEW(nodempc,ITG,3*memmpc_);
744  RENEW(coefmpc,double,memmpc_);
745  for(k=index;k<memmpc_;k++){
746  nodempc[3*k-1]=k+1;
747  }
748  nodempc[3*memmpc_-1]=0;
749 
750  /* creating the A.n MPC's */
751 
752  FORTRAN(createinterfacempcs,(imastnode,xmastnor,&nmastnode,
753  ikmpc,ilmpc,nmpc,ipompc,nodempc,coefmpc,labmpc,&mpcfree,
754  ikboun,nboun));
755 
756  SFREE(imastnode);SFREE(xmastnor);
757  }
758 
759  /* determining the new matrix structure */
760 
761  remastructem(ipompc,&coefmpc,&nodempc,nmpc,
762  &mpcfree,nodeboun,ndirboun,nboun,ikmpc,ilmpc,ikboun,ilboun,
763  labmpc,nk,&memmpc_,&icascade,&maxlenmpc,
764  kon,ipkon,lakon,ne,nactdof,icol,jq,&irow,isolver,
765  neq,nzs,nmethod,&f,&fext,&b,&aux2,&fini,&fextini,
766  &adb,&aub,ithermal,iperturb,mass,mi,ielmat,elcon,
767  ncmat_,ntmat_,inomat,network);
768 
769 /**************************************************************/
770 /* starting the loop over the increments */
771 /**************************************************************/
772 
773  /* saving the distributed loads (volume heating will be
774  added because of Joule heating) */
775 
776  if(*ithermal==3){
777  nloadref=*nload;
778  NNEW(nelemloadref,ITG,2**nload);
779  if(*nam>0) NNEW(iamloadref,ITG,2**nload);
780  NNEW(sideloadref,char,20**nload);
781 
782  memcpy(&nelemloadref[0],&nelemload[0],sizeof(ITG)*2**nload);
783  if(*nam>0) memcpy(&iamloadref[0],&iamload[0],sizeof(ITG)*2**nload);
784  memcpy(&sideloadref[0],&sideload[0],sizeof(char)*20**nload);
785 
786  /* generating new fields; ne2 is the number of elements
787  in domain 2 = A,V-domain (the only domain with heating) */
788 
789  (*nload_)+=ne2;
790  RENEW(nelemload,ITG,2**nload_);
791  if(*nam>0) RENEW(iamload,ITG,2**nload_);
792  RENEW(xloadact,double,2**nload_);
793  RENEW(sideload,char,20**nload_);
794  }
795 
796  if((*ithermal==1)||(*ithermal>=3)){
797  NNEW(t1ini,double,*nk);
798  NNEW(t1act,double,*nk);
799  for(k=0;k<*nk;++k){t1act[k]=t1old[k];}
800  }
801 
802  newstep=1;
803 
804  while(1.-theta>1.e-6){
805 
806  if(icutb==0){
807 
808  /* previous increment converged: update the initial values */
809 
810  iinc++;
811  jprint++;
812 
813  /* vold is copied into vini */
814 
815  memcpy(&vini[0],&vold[0],sizeof(double)*mt**nk);
816 
817  for(k=0;k<*nboun;++k){xbounini[k]=xbounact[k];}
818  if((*ithermal==1)||(*ithermal>=3)){
819  for(k=0;k<*nk;++k){t1ini[k]=t1act[k];}
820  }
821  for(k=0;k<neq[1];++k){
822  fini[k]=f[k];
823  }
824  if(*nmethod==4){
825  for(k=0;k<mt**nk;++k){
826  veini[k]=veold[k];
827  }
828  for(k=0;k<neq[1];++k){
829  fextini[k]=fext[k];
830  }
831  }
832  }
833 
834  /* check for max. # of increments */
835 
836  if(iinc>*jmax){
837  printf(" *ERROR: max. # of increments reached\n\n");
838  FORTRAN(stop,());
839  }
840  printf(" increment %" ITGFORMAT " attempt %" ITGFORMAT " \n",iinc,icutb+1);
841  printf(" increment size= %e\n",dtheta**tper);
842  printf(" sum of previous increments=%e\n",theta**tper);
843  printf(" actual step time=%e\n",(theta+dtheta)**tper);
844  printf(" actual total time=%e\n\n",*ttime+(theta+dtheta)**tper);
845 
846  printf(" iteration 1\n\n");
847 
848  qamold[0]=qam[0];
849  qamold[1]=qam[1];
850 
851  /* determining the actual loads at the end of the new increment*/
852 
853  reltime=theta+dtheta;
854  time=reltime**tper;
855  dtime=dtheta**tper;
856 
857  /* restoring the distributed loading before adding the
858  Joule heating */
859 
860  if(*ithermal==3){
861  *nload=nloadref;
862  DMEMSET(nelemload,0,2**nload_,0);
863  memcpy(&nelemload[0],&nelemloadref[0],sizeof(ITG)*2**nload);
864  if(*nam>0){
865  DMEMSET(iamload,0,2**nload_,0);
866  memcpy(&iamload[0],&iamloadref[0],sizeof(ITG)*2**nload);
867  }
868  DMEMSET(xloadact,0,2**nload_,0.);
869  DMEMSET(sideload,0,'\0',0.);memcpy(&sideload[0],&sideloadref[0],sizeof(char)*20**nload);
870  }
871 
872  /* determining the actual loading */
873 
874 // for(i=0;i<3**nk;i++){h0[i]/=h0scale;}
875  FORTRAN(tempload_em,(xforcold,xforc,xforcact,iamforc,nforc,xloadold,xload,
876  xloadact,iamload,nload,ibody,xbody,nbody,xbodyold,xbodyact,
877  t1old,t1,t1act,iamt1,nk,amta,
878  namta,nam,ampli,&time,&reltime,ttime,&dtime,ithermal,nmethod,
879  xbounold,xboun,xbounact,iamboun,nboun,
880  nodeboun,ndirboun,nodeforc,ndirforc,istep,&iinc,
881  co,vold,itg,&ntg,amname,ikboun,ilboun,nelemload,sideload,mi,
882  ntrans,trab,inotr,veold,integerglob,doubleglob,tieset,istartset,
883  iendset,ialset,ntie,nmpc,ipompc,ikmpc,ilmpc,nodempc,coefmpc,
884  &h0scale,inomat,ipobody,iponoel,inoel));
885  for(i=0;i<3**nk;i++){h0[i]=h0ref[i]*h0scale;}
886 
887  for(i=0;i<3;i++){cam[i]=0.;}for(i=3;i<5;i++){cam[i]=0.5;}
888  if(*ithermal>1){radflowload(itg,ieg,&ntg,&ntr,adrad,aurad,bcr,ipivr,
889  ac,bc,nload,sideload,nelemload,xloadact,lakon,ipiv,ntmat_,vold,
890  shcon,nshcon,ipkon,kon,co,
891  kontri,&ntri,nloadtr,tarea,tenv,physcon,erad,&adview,&auview,
892  nflow,ikboun,xbounact,nboun,ithermal,&iinc,&iit,
893  cs,mcs,inocs,&ntrit,nk,fenv,istep,&dtime,ttime,&time,ilboun,
894  ikforc,ilforc,xforcact,nforc,cam,ielmat,&nteq,prop,ielprop,
895  nactdog,nacteq,nodeboun,ndirboun,network,
896  rhcon,nrhcon,ipobody,ibody,xbodyact,nbody,iviewfile,jobnamef,
897  ctrl,xloadold,&reltime,nmethod,set,mi,istartset,iendset,ialset,nset,
898  ineighe,nmpc,nodempc,ipompc,coefmpc,labmpc,&iemchange,nam,iamload,
899  jqrad,irowrad,&nzsrad,icolrad,ne,iaxial,qa,cocon,ncocon,iponoel,
900  inoel,nprop,amname,namta,amta);
901  }
902 
903  /* prediction of the next solution (only for temperature) */
904 
905  NNEW(v,double,mt**nk);
906 
907 // if(*ithermal>2){
908  prediction_em(uam,nmethod,&bet,&gam,&dtime,ithermal,nk,veold,v,
909  &iinc,&idiscon,vold,nactdof,mi);
910 // }
911 
912  NNEW(fn,double,mt**nk);
913 
914  iout=-1;
915  iperturb_sav[0]=iperturb[0];
916  iperturb_sav[1]=iperturb[1];
917 
918  /* first iteration in first increment: heat tangent */
919 
920  NNEW(inum,ITG,*nk);
921  resultsinduction(co,nk,kon,ipkon,lakon,ne,v,stn,inum,
922  elcon,nelcon,rhcon,nrhcon,alcon,nalcon,alzero,ielmat,
923  ielorien,norien,orab,ntmat_,t0,t1act,ithermal,
924  prestr,iprestr,filab,eme,emn,een,iperturb,
925  f,fn,nactdof,&iout,qa,vold,b,nodeboun,
926  ndirboun,xbounact,nboun,ipompc,
927  nodempc,coefmpc,labmpc,nmpc,nmethod,cam,&neq[1],veold,accold,
928  &bet,&gam,&dtime,&time,ttime,plicon,nplicon,plkcon,nplkcon,
929  xstateini,xstiff,xstate,npmat_,epn,matname,mi,&ielas,
930  &icmd,ncmat_,nstate_,sti,vini,ikboun,ilboun,ener,enern,
931  emeini,xstaten,eei,enerini,cocon,ncocon,set,nset,istartset,
932  iendset,ialset,nprint,prlab,prset,qfx,qfn,trab,inotr,ntrans,
933  fmpc,nelemload,nload,ikmpc,ilmpc,istep,&iinc,springarea,
934  &reltime,ne,xforc,nforc,thicke,shcon,nshcon,
935  sideload,xloadact,xloadold,&icfd,inomat,h0,islavnode,
936  nslavnode,ntie,ielprop,prop,iactive,energyini,energy,
937  iponoel,inoel,orname,network,ipobody,xbodyact,ibody);
938  SFREE(inum);
939 
940  /* the calculation of the electromagnetic fields is (quasi)linear,
941  i.e. the solution of the equations is the fields;
942  only the temperature calculation is nonlinear,
943  i.e. the solution of the equations is a differential temperature */
944 
945  memcpy(&vold[0],&v[0],sizeof(double)*mt**nk);
946 
947  iout=0;
948 
949  SFREE(fn);SFREE(v);
950 
951  /***************************************************************/
952  /* iteration counter and start of the loop over the iterations */
953  /***************************************************************/
954 
955  iit=1;
956  icntrl=0;
957  ctrl[0]=i0ref;ctrl[1]=irref;ctrl[3]=icref;
958 
959  while(icntrl==0){
960 
961  if(iit!=1){
962 
963  printf(" iteration %" ITGFORMAT "\n\n",iit);
964 
965  /* restoring the distributed loading before adding the
966  Joule heating */
967 
968  if(*ithermal==3){
969  *nload=nloadref;
970  DMEMSET(nelemload,0,2**nload_,0);
971  memcpy(&nelemload[0],&nelemloadref[0],sizeof(ITG)*2**nload);
972  if(*nam>0){
973  DMEMSET(iamload,0,2**nload_,0);
974  memcpy(&iamload[0],&iamloadref[0],sizeof(ITG)*2**nload);
975  }
976  DMEMSET(xloadact,0,2**nload_,0.);
977  DMEMSET(sideload,0,20**nload_,'\0');memcpy(&sideload[0],&sideloadref[0],sizeof(char)*20**nload);
978  }
979 
980  /* determining the actual loading */
981 
982 // for(i=0;i<3**nk;i++){h0[i]/=h0scale;}
983  FORTRAN(tempload_em,(xforcold,xforc,xforcact,iamforc,nforc,
984  xloadold,xload,
985  xloadact,iamload,nload,ibody,xbody,nbody,xbodyold,xbodyact,
986  t1old,t1,t1act,iamt1,nk,amta,
987  namta,nam,ampli,&time,&reltime,ttime,&dtime,ithermal,nmethod,
988  xbounold,xboun,xbounact,iamboun,nboun,
989  nodeboun,ndirboun,nodeforc,ndirforc,istep,&iinc,
990  co,vold,itg,&ntg,amname,ikboun,ilboun,nelemload,sideload,mi,
991  ntrans,trab,inotr,veold,integerglob,doubleglob,tieset,istartset,
992  iendset,ialset,ntie,nmpc,ipompc,ikmpc,ilmpc,nodempc,coefmpc,
993  &h0scale,inomat,ipobody,iponoel,inoel));
994  for(i=0;i<3**nk;i++){h0[i]=h0ref[i]*h0scale;}
995 
996  for(i=0;i<3;i++){cam[i]=0.;}for(i=3;i<5;i++){cam[i]=0.5;}
997  if(*ithermal>1){radflowload(itg,ieg,&ntg,&ntr,adrad,aurad,
998  bcr,ipivr,ac,bc,nload,sideload,nelemload,xloadact,lakon,ipiv,
999  ntmat_,vold,shcon,nshcon,ipkon,kon,co,
1000  kontri,&ntri,nloadtr,tarea,tenv,physcon,erad,&adview,&auview,
1001  nflow,ikboun,xbounact,nboun,ithermal,&iinc,&iit,
1002  cs,mcs,inocs,&ntrit,nk,fenv,istep,&dtime,ttime,&time,ilboun,
1003  ikforc,ilforc,xforcact,nforc,cam,ielmat,&nteq,prop,ielprop,
1004  nactdog,nacteq,nodeboun,ndirboun,network,
1005  rhcon,nrhcon,ipobody,ibody,xbodyact,nbody,iviewfile,jobnamef,
1006  ctrl,xloadold,&reltime,nmethod,set,mi,istartset,iendset,ialset,
1007  nset,ineighe,nmpc,nodempc,ipompc,coefmpc,labmpc,&iemchange,nam,
1008  iamload,jqrad,irowrad,&nzsrad,icolrad,ne,iaxial,qa,cocon,
1009  ncocon,iponoel,inoel,nprop,amname,namta,amta);
1010  }
1011 
1012  }
1013 
1014  /* add Joule heating */
1015 
1016  if(*ithermal==3){
1017  NNEW(idefload,ITG,*nload_);
1018  DMEMSET(idefload,0,*nload_,1);
1019  FORTRAN(jouleheating,(ipkon,lakon,kon,co,elcon,nelcon,
1020  mi,ne,sti,ielmat,nelemload,sideload,xloadact,nload,nload_,
1021  iamload,nam,idefload,ncmat_,ntmat_,
1022  alcon,nalcon,ithermal,vold,t1));
1023  SFREE(idefload);
1024  }
1025 
1026  if(*ithermal==3){
1027  for(k=0;k<*nk;++k){t1act[k]=vold[mt*k];}
1028  }
1029 
1030  /* calculating the local stiffness matrix and external loading */
1031 
1032  NNEW(ad,double,neq[1]);
1033  NNEW(au,double,nzs[1]);
1034 
1035  FORTRAN(mafillem,(co,nk,kon,ipkon,lakon,ne,nodeboun,ndirboun,
1036  xbounact,nboun,
1037  ipompc,nodempc,coefmpc,nmpc,nodeforc,ndirforc,xforcact,
1038  nforc,nelemload,sideload,xloadact,nload,xbodyact,ipobody,
1039  nbody,cgr,ad,au,fext,nactdof,icol,jq,irow,neq,nzl,
1040  nmethod,ikmpc,ilmpc,ikboun,ilboun,
1041  elcon,nelcon,rhcon,nrhcon,alcon,nalcon,alzero,
1042  ielmat,ielorien,norien,orab,ntmat_,
1043  t0,t1act,ithermal,prestr,iprestr,vold,iperturb,sti,
1044  nzs,stx,adb,aub,iexpl,plicon,nplicon,plkcon,nplkcon,
1045  xstiff,npmat_,&dtime,matname,mi,
1046  ncmat_,mass,&stiffness,&buckling,&rhsi,&intscheme,
1047  physcon,shcon,nshcon,cocon,ncocon,ttime,&time,istep,&iinc,
1048  &coriolis,ibody,xloadold,&reltime,veold,springarea,nstate_,
1049  xstateini,xstate,thicke,integerglob,doubleglob,
1050  tieset,istartset,iendset,ialset,ntie,&nasym,iactive,h0,
1051  pslavsurf,pmastsurf,&mortar,clearini,ielprop,prop,
1052  iponoel,inoel,network));
1053 
1054  iperturb[0]=iperturb_sav[0];
1055  iperturb[1]=iperturb_sav[1];
1056 
1057  /* calculating the residual (f is only for the temperature
1058  nonzero) */
1059 
1060  calcresidual_em(nmethod,neq,b,fext,f,iexpl,nactdof,aux1,aux2,vold,
1061  vini,&dtime,accold,nk,adb,aub,jq,irow,nzl,alpha,fextini,fini,
1062  islavnode,nslavnode,&mortar,ntie,f_cm,f_cs,mi,
1063  nzs,&nasym,ithermal);
1064 
1065  newstep=0;
1066 
1067  if(*nmethod==0){
1068 
1069  /* error occurred in mafill: storing the geometry in frd format */
1070 
1071  *nmethod=0;
1072  ++*kode;
1073  NNEW(inum,ITG,*nk);for(k=0;k<*nk;k++) inum[k]=1;
1074 
1075  ptime=*ttime+time;
1076  frd(co,nk,kon,ipkon,lakon,ne,v,stn,inum,nmethod,
1077  kode,filab,een,t1,fn,&ptime,epn,ielmat,matname,enern,xstaten,
1078  nstate_,istep,&iinc,ithermal,qfn,&mode,&noddiam,trab,inotr,
1079  ntrans,orab,ielorien,norien,description,ipneigh,neigh,
1080  mi,sti,vr,vi,stnr,stni,vmax,stnmax,&ngraph,veold,ener,ne,
1081  cs,set,nset,istartset,iendset,ialset,eenmax,fnr,fni,emn,
1082  thicke,jobnamec,output,qfx,cdn,&mortar,cdnr,cdni,nmat);
1083 
1084  }
1085 
1086  /* implicit step (static or dynamic */
1087 
1088  if(*nmethod==4){
1089 
1090  /* electromagnetic part */
1091 
1092  if(*ithermal!=2){
1093  for(k=0;k<neq[0];++k){
1094  ad[k]=adb[k]/dtime+ad[k];
1095  }
1096  for(k=0;k<nzs[0];++k){
1097  au[k]=aub[k]/dtime+au[k];
1098  }
1099 
1100  /* upper triangle of asymmetric matrix */
1101 
1102  if(nasym>0){
1103  for(k=nzs[2];k<nzs[2]+nzs[0];++k){
1104  au[k]=aub[k]/dtime+au[k];
1105  }
1106  }
1107  }
1108 
1109  /* thermal part */
1110 
1111  if(*ithermal>1){
1112  for(k=neq[0];k<neq[1];++k){
1113  ad[k]=adb[k]/dtime+ad[k];
1114  }
1115  for(k=nzs[0];k<nzs[1];++k){
1116  au[k]=aub[k]/dtime+au[k];
1117  }
1118 
1119  /* upper triangle of asymmetric matrix */
1120 
1121  if(nasym>0){
1122  for(k=nzs[2]+nzs[0];k<nzs[2]+nzs[1];++k){
1123  au[k]=aub[k]/dtime+au[k];
1124  }
1125  }
1126  }
1127  }
1128 
1129  NNEW(adaux,double,neq[2]);
1130  FORTRAN(preconditioning,(ad,au,b,&neq[1],irow,jq,adaux));
1131 
1132 
1133  if(*isolver==0){
1134 #ifdef SPOOLES
1135  if(*ithermal<2){
1136  spooles(ad,au,adb,aub,&sigma,b,icol,irow,&neq[0],&nzs[0],
1137  &symmetryflag,&inputformat,&nzs[2]);
1138  }else{
1139  spooles(ad,au,adb,aub,&sigma,b,icol,irow,&neq[1],&nzs[1],
1140  &symmetryflag,&inputformat,&nzs[2]);
1141  }
1142 #else
1143  printf(" *ERROR in nonlingeo: the SPOOLES library is not linked\n\n");
1144  FORTRAN(stop,());
1145 #endif
1146  }
1147  else if((*isolver==2)||(*isolver==3)){
1148  if(nasym>0){
1149  printf(" *ERROR in nonlingeo: the iterative solver cannot be used for asymmetric matrices\n\n");
1150  FORTRAN(stop,());
1151  }
1152  preiter(ad,&au,b,&icol,&irow,&neq[1],&nzs[1],isolver,iperturb);
1153  }
1154  else if(*isolver==4){
1155 #ifdef SGI
1156  if(nasym>0){
1157  printf(" *ERROR in nonlingeo: the SGI solver cannot be used for asymmetric matrices\n\n");
1158  FORTRAN(stop,());
1159  }
1160  token=1;
1161  if(*ithermal<2){
1162  sgi_main(ad,au,adb,aub,&sigma,b,icol,irow,&neq[0],&nzs[0],token);
1163  }else{
1164  sgi_main(ad,au,adb,aub,&sigma,b,icol,irow,&neq[1],&nzs[1],token);
1165  }
1166 #else
1167  printf(" *ERROR in nonlingeo: the SGI library is not linked\n\n");
1168  FORTRAN(stop,());
1169 #endif
1170  }
1171  else if(*isolver==5){
1172  if(nasym>0){
1173  printf(" *ERROR in nonlingeo: the TAUCS solver cannot be used for asymmetric matrices\n\n");
1174  FORTRAN(stop,());
1175  }
1176 #ifdef TAUCS
1177  tau(ad,&au,adb,aub,&sigma,b,icol,&irow,&neq[1],&nzs[1]);
1178 #else
1179  printf(" *ERROR in nonlingeo: the TAUCS library is not linked\n\n");
1180  FORTRAN(stop,());
1181 #endif
1182  }
1183  else if(*isolver==7){
1184 #ifdef PARDISO
1185  if(*ithermal<2){
1186  pardiso_main(ad,au,adb,aub,&sigma,b,icol,irow,&neq[0],&nzs[0],
1187  &symmetryflag,&inputformat,jq,&nzs[2]);
1188  }else{
1189  pardiso_main(ad,au,adb,aub,&sigma,b,icol,irow,&neq[1],&nzs[1],
1190  &symmetryflag,&inputformat,jq,&nzs[2]);
1191  }
1192 #else
1193  printf(" *ERROR in nonlingeo: the PARDISO library is not linked\n\n");
1194  FORTRAN(stop,());
1195 #endif
1196  }
1197 
1198  for(i=0;i<neq[1];i++){b[i]*=adaux[i];}
1199  SFREE(adaux);
1200 
1201  /* calculating the electromagnetic fields and temperatures
1202  only the temperature calculation is differential */
1203 
1204  NNEW(v,double,mt**nk);
1205  memcpy(&v[0],&vold[0],sizeof(double)*mt**nk);
1206 
1207  NNEW(fn,double,mt**nk);
1208 
1209  NNEW(inum,ITG,*nk);
1210  resultsinduction(co,nk,kon,ipkon,lakon,ne,v,stn,inum,
1211  elcon,nelcon,rhcon,nrhcon,alcon,nalcon,alzero,ielmat,
1212  ielorien,norien,orab,ntmat_,t0,t1act,ithermal,
1213  prestr,iprestr,filab,eme,emn,een,iperturb,
1214  f,fn,nactdof,&iout,qa,vold,b,nodeboun,
1215  ndirboun,xbounact,nboun,ipompc,
1216  nodempc,coefmpc,labmpc,nmpc,nmethod,cam,&neq[1],veold,accold,
1217  &bet,&gam,&dtime,&time,ttime,plicon,nplicon,plkcon,nplkcon,
1218  xstateini,xstiff,xstate,npmat_,epn,matname,mi,&ielas,
1219  &icmd,ncmat_,nstate_,sti,vini,ikboun,ilboun,ener,enern,
1220  emeini,xstaten,eei,enerini,cocon,ncocon,set,nset,istartset,
1221  iendset,ialset,nprint,prlab,prset,qfx,qfn,trab,inotr,ntrans,
1222  fmpc,nelemload,nload,ikmpc,ilmpc,istep,&iinc,springarea,
1223  &reltime,ne,xforc,nforc,thicke,shcon,nshcon,
1224  sideload,xloadact,xloadold,&icfd,inomat,h0,islavnode,
1225  nslavnode,ntie,ielprop,prop,iactive,energyini,energy,
1226  iponoel,inoel,orname,network,ipobody,xbodyact,ibody);
1227  SFREE(inum);
1228 
1229  SFREE(ad);SFREE(au);
1230 
1231  if(*ithermal!=2){
1232  if(cam[0]>uam[0]){uam[0]=cam[0];}
1233  if(qau<1.e-10){
1234  if(qa[0]>ea*qam[0]){qam[0]=(qamold[0]*jnz+qa[0])/(jnz+1);}
1235  else {qam[0]=qamold[0];}
1236  }
1237  }
1238  if(*ithermal>1){
1239  if(cam[1]>uam[1]){uam[1]=cam[1];}
1240  if(qau<1.e-10){
1241  if(qa[1]>ea*qam[1]){qam[1]=(qamold[1]*jnz+qa[1])/(jnz+1);}
1242  else {qam[1]=qamold[1];}
1243  }
1244  }
1245 
1246  memcpy(&vold[0],&v[0],sizeof(double)*mt**nk);
1247 
1248  SFREE(v);SFREE(fn);
1249 
1250  /* calculating the residual */
1251 
1252  calcresidual_em(nmethod,neq,b,fext,f,iexpl,nactdof,aux1,aux2,vold,
1253  vini,&dtime,accold,nk,adb,aub,jq,irow,nzl,alpha,fextini,fini,
1254  islavnode,nslavnode,&mortar,ntie,f_cm,f_cs,mi,
1255  nzs,&nasym,ithermal);
1256 
1257  /* calculating the maximum residual (only thermal part)*/
1258 
1259  for(k=0;k<2;++k){
1260  ram2[k]=ram1[k];
1261  ram1[k]=ram[k];
1262  ram[k]=0.;
1263  }
1264 
1265  if(*ithermal!=2) ram[0]=0.;
1266 
1267  if(*ithermal>1){
1268  for(k=neq[0];k<neq[1];++k){
1269  err=fabs(b[k]);
1270  if(err>ram[1]){ram[1]=err;ram[3]=k+0.5;}
1271  }
1272  }
1273 
1274  /* printing residuals */
1275 
1276  if(*ithermal>1){
1277  if(ram[1]<1.e-6) ram[1]=0.;
1278  printf(" average flux= %f\n",qa[1]);
1279  printf(" time avg. flux= %f\n",qam[1]);
1280  if((ITG)((double)nactdofinv[(ITG)ram[3]]/mt)+1==0){
1281  printf(" largest residual flux= %f\n",
1282  ram[1]);
1283  }else{
1284  inode=(ITG)((double)nactdofinv[(ITG)ram[3]]/mt)+1;
1285  idir=nactdofinv[(ITG)ram[3]]-mt*(inode-1);
1286  printf(" largest residual flux= %f in node %" ITGFORMAT " and dof %" ITGFORMAT "\n",ram[1],inode,idir);
1287  }
1288  printf(" largest increment of temp= %e\n",uam[1]);
1289  if((ITG)cam[4]==0){
1290  printf(" largest correction to temp= %e\n\n",
1291  cam[1]);
1292  }else{
1293  inode=(ITG)((double)nactdofinv[(ITG)cam[4]]/mt)+1;
1294  idir=nactdofinv[(ITG)cam[4]]-mt*(inode-1);
1295  printf(" largest correction to temp= %e in node %" ITGFORMAT " and dof %" ITGFORMAT "\n\n",cam[1],inode,idir);
1296  }
1297  }
1298 
1299  /* athermal electromagnetic calculations are linear:
1300  set iit=2 to force convergence */
1301 
1302  if(*ithermal<=1) iit=2;
1303 
1304  // MPADD: need for fake energy values!
1305  double energy[4] = {0, 0, 0, 0};
1306  double allwk = 0.0;
1307  double energyref = 0.0;
1308  double emax, enres,enetoll, reswk, dampwk, allwkini;
1309 
1310  neini=*ne;
1311  checkconvergence(co,nk,kon,ipkon,lakon,ne,stn,nmethod,
1312  kode,filab,een,t1act,&time,epn,ielmat,matname,enern,
1313  xstaten,nstate_,istep,&iinc,iperturb,ener,mi,output,
1314  ithermal,qfn,&mode,&noddiam,trab,inotr,ntrans,orab,
1315  ielorien,norien,description,sti,&icutb,&iit,&dtime,qa,
1316  vold,qam,ram1,ram2,ram,cam,uam,&ntg,ttime,&icntrl,
1317  &theta,&dtheta,veold,vini,idrct,tper,&istab,tmax,
1318  nactdof,b,tmin,ctrl,amta,namta,itpamp,&inext,&dthetaref,
1319  &itp,&jprint,jout,&uncoupled,t1,&iitterm,nelemload,
1320  nload,nodeboun,nboun,itg,ndirboun,&deltmx,&iflagact,
1321  set,nset,istartset,iendset,ialset,emn,thicke,jobnamec,
1322  &mortar,nmat,ielprop,prop,&ialeatoric,&kscale,
1323  energy, &allwk, &energyref,&emax, &enres, &enetoll, //MPADD
1324  energyini, &allwkini ,&allwk, &reswk, &ne0, &ne0, &dampwk, //MPADD
1325  &dampwk, energy); //MPADD
1326  }
1327 
1328  /*********************************************************/
1329  /* end of the iteration loop */
1330  /*********************************************************/
1331 
1332  /* icutb=0 means that the iterations in the increment converged,
1333  icutb!=0 indicates that the increment has to be reiterated with
1334  another increment size (dtheta) */
1335 
1336  if(((qa[0]>ea*qam[0])||(qa[1]>ea*qam[1]))&&(icutb==0)){jnz++;}
1337  iit=0;
1338 
1339  if(icutb!=0){
1340  memcpy(&vold[0],&vini[0],sizeof(double)*mt**nk);
1341 
1342  for(k=0;k<*nboun;++k){xbounact[k]=xbounini[k];}
1343  if((*ithermal==1)||(*ithermal>=3)){
1344  for(k=0;k<*nk;++k){t1act[k]=t1ini[k];}
1345  }
1346  for(k=0;k<neq[1];++k){
1347  f[k]=fini[k];
1348  }
1349  if(*nmethod==4){
1350  for(k=0;k<mt**nk;++k){
1351  veold[k]=veini[k];
1352  }
1353  for(k=0;k<neq[1];++k){
1354  fext[k]=fextini[k];
1355  }
1356  }
1357 
1358  qam[0]=qamold[0];
1359  qam[1]=qamold[1];
1360  }
1361 
1362  if((jout[0]==jprint)&&(icutb==0)){
1363 
1364  jprint=0;
1365 
1366  /* calculating the displacements and the stresses and storing */
1367  /* the results in frd format */
1368 
1369  NNEW(v,double,mt**nk);
1370  NNEW(fn,double,mt**nk);
1371  if(*ithermal>1) NNEW(qfn,double,3**nk);
1372  if((strcmp1(&filab[3741],"EMFE")==0)||
1373  (strcmp1(&filab[3828],"EMFB")==0)) NNEW(stn,double,6**nk);
1374  NNEW(inum,ITG,*nk);
1375 
1376  memcpy(&v[0],&vold[0],sizeof(double)*mt**nk);
1377 
1378  iout=2;
1379  icmd=3;
1380 
1381  resultsinduction(co,nk,kon,ipkon,lakon,ne,v,stn,inum,
1382  elcon,nelcon,rhcon,nrhcon,alcon,nalcon,alzero,ielmat,
1383  ielorien,norien,orab,ntmat_,t0,t1act,ithermal,
1384  prestr,iprestr,filab,eme,emn,een,iperturb,
1385  f,fn,nactdof,&iout,qa,vold,b,nodeboun,
1386  ndirboun,xbounact,nboun,ipompc,
1387  nodempc,coefmpc,labmpc,nmpc,nmethod,cam,&neq[1],veold,accold,
1388  &bet,&gam,&dtime,&time,ttime,plicon,nplicon,plkcon,nplkcon,
1389  xstateini,xstiff,xstate,npmat_,epn,matname,mi,&ielas,&icmd,
1390  ncmat_,nstate_,sti,vini,ikboun,ilboun,ener,enern,emeini,
1391  xstaten,eei,enerini,cocon,ncocon,set,nset,istartset,iendset,
1392  ialset,nprint,prlab,prset,qfx,qfn,trab,inotr,ntrans,fmpc,
1393  nelemload,nload,ikmpc,ilmpc,istep,&iinc,springarea,
1394  &reltime,ne,xforc,nforc,thicke,shcon,nshcon,
1395  sideload,xloadact,xloadold,&icfd,inomat,h0,islavnode,
1396  nslavnode,ntie,ielprop,prop,iactive,energyini,energy,
1397  iponoel,inoel,orname,network,ipobody,xbodyact,ibody);
1398 
1399  memcpy(&vold[0],&v[0],sizeof(double)*mt**nk);
1400 
1401  iout=0;
1402  icmd=0;
1403 // FORTRAN(networkinum,(ipkon,inum,kon,lakon,ne,itg,&ntg));
1404 
1405  ++*kode;
1406  if(*mcs!=0){
1407  ptime=*ttime+time;
1408  frdcyc(co,nk,kon,ipkon,lakon,ne,v,stn,inum,nmethod,kode,filab,een,
1409  t1act,fn,&ptime,epn,ielmat,matname,cs,mcs,nkon,enern,xstaten,
1410  nstate_,istep,&iinc,iperturb,ener,mi,output,ithermal,qfn,
1411  ialset,istartset,iendset,trab,inotr,ntrans,orab,ielorien,
1412  norien,stx,veold,&noddiam,set,nset,emn,thicke,jobnamec,ne,
1413  cdn,&mortar,nmat,qfx);
1414  }else{
1415 
1416  ptime=*ttime+time;
1417  frd(co,nk,kon,ipkon,lakon,ne,v,stn,inum,nmethod,
1418  kode,filab,een,t1act,fn,&ptime,epn,ielmat,matname,
1419  enern,xstaten,
1420  nstate_,istep,&iinc,ithermal,qfn,&mode,&noddiam,trab,inotr,
1421  ntrans,orab,ielorien,norien,description,ipneigh,neigh,
1422  mi,stx,vr,vi,stnr,stni,vmax,stnmax,&ngraph,veold,ener,ne,
1423  cs,set,nset,istartset,iendset,ialset,eenmax,fnr,fni,emn,
1424  thicke,jobnamec,output,qfx,cdn,&mortar,cdnr,cdni,nmat);
1425 
1426  }
1427 
1428  SFREE(v);SFREE(fn);SFREE(inum);
1429  if(*ithermal>1){SFREE(qfn);}
1430  if((strcmp1(&filab[3741],"EMFE")==0)||
1431  (strcmp1(&filab[3828],"EMFB")==0)) SFREE(stn);
1432 
1433  }
1434 
1435  }
1436 
1437  /*********************************************************/
1438  /* end of the increment loop */
1439  /*********************************************************/
1440 
1441  if(jprint!=0){
1442 
1443  /* calculating the displacements and the stresses and storing
1444  the results in frd format */
1445 
1446  NNEW(v,double,mt**nk);
1447  NNEW(fn,double,mt**nk);
1448  if(*ithermal>1) NNEW(qfn,double,3**nk);
1449  if(strcmp1(&filab[3741],"EMF ")==0) NNEW(stn,double,6**nk);
1450  NNEW(inum,ITG,*nk);
1451 
1452  memcpy(&v[0],&vold[0],sizeof(double)*mt**nk);
1453  iout=2;
1454  icmd=3;
1455 
1456  resultsinduction(co,nk,kon,ipkon,lakon,ne,v,stn,inum,
1457  elcon,nelcon,rhcon,nrhcon,alcon,nalcon,alzero,ielmat,
1458  ielorien,norien,orab,ntmat_,t0,t1act,ithermal,
1459  prestr,iprestr,filab,eme,emn,een,iperturb,
1460  f,fn,nactdof,&iout,qa,vold,b,nodeboun,
1461  ndirboun,xbounact,nboun,ipompc,
1462  nodempc,coefmpc,labmpc,nmpc,nmethod,cam,&neq[1],veold,accold,
1463  &bet,&gam,&dtime,&time,ttime,plicon,nplicon,plkcon,nplkcon,
1464  xstateini,xstiff,xstate,npmat_,epn,matname,mi,&ielas,&icmd,
1465  ncmat_,nstate_,sti,vini,ikboun,ilboun,ener,enern,emeini,
1466  xstaten,eei,enerini,cocon,ncocon,set,nset,istartset,iendset,
1467  ialset,nprint,prlab,prset,qfx,qfn,trab,inotr,ntrans,fmpc,
1468  nelemload,nload,ikmpc,ilmpc,istep,&iinc,springarea,
1469  &reltime,ne,xforc,nforc,thicke,shcon,nshcon,
1470  sideload,xloadact,xloadold,&icfd,inomat,h0,islavnode,
1471  nslavnode,ntie,ielprop,prop,iactive,energyini,energy,
1472  iponoel,inoel,orname,network,ipobody,xbodyact,ibody);
1473 
1474  memcpy(&vold[0],&v[0],sizeof(double)*mt**nk);
1475 
1476  iout=0;
1477  icmd=0;
1478 // FORTRAN(networkinum,(ipkon,inum,kon,lakon,ne,itg,&ntg));
1479 
1480  ++*kode;
1481  if(*mcs>0){
1482  ptime=*ttime+time;
1483  frdcyc(co,nk,kon,ipkon,lakon,ne,v,stn,inum,nmethod,kode,filab,een,
1484  t1act,fn,&ptime,epn,ielmat,matname,cs,mcs,nkon,enern,xstaten,
1485  nstate_,istep,&iinc,iperturb,ener,mi,output,ithermal,qfn,
1486  ialset,istartset,iendset,trab,inotr,ntrans,orab,ielorien,
1487  norien,stx,veold,&noddiam,set,nset,emn,thicke,jobnamec,ne,
1488  cdn,&mortar,nmat,qfx);
1489  }else{
1490 
1491  ptime=*ttime+time;
1492  frd(co,nk,kon,ipkon,lakon,ne,v,stn,inum,nmethod,
1493  kode,filab,een,t1act,fn,&ptime,epn,ielmat,matname,enern,xstaten,
1494  nstate_,istep,&iinc,ithermal,qfn,&mode,&noddiam,trab,inotr,
1495  ntrans,orab,ielorien,norien,description,ipneigh,neigh,
1496  mi,stx,vr,vi,stnr,stni,vmax,stnmax,&ngraph,veold,ener,ne,
1497  cs,set,nset,istartset,iendset,ialset,eenmax,fnr,fni,emn,
1498  thicke,jobnamec,output,qfx,cdn,&mortar,cdnr,cdni,nmat);
1499 
1500  }
1501 
1502  SFREE(v);SFREE(fn);SFREE(inum);
1503  if(*ithermal>1){SFREE(qfn);}
1504  if(strcmp1(&filab[3741],"EMF ")==0) SFREE(stn);
1505 
1506  }
1507 
1508  /* restoring the distributed loading */
1509 
1510  if(*ithermal==3){
1511  *nload=nloadref;
1512  (*nload_)-=ne2;
1513  RENEW(nelemload,ITG,2**nload);memcpy(&nelemload[0],&nelemloadref[0],sizeof(ITG)*2**nload);
1514  if(*nam>0){
1515  RENEW(iamload,ITG,2**nload);
1516  memcpy(&iamload[0],&iamloadref[0],sizeof(ITG)*2**nload);
1517  }
1518  RENEW(sideload,char,20**nload);memcpy(&sideload[0],&sideloadref[0],sizeof(char)*20**nload);
1519 
1520  /* freeing the temporary fields */
1521 
1522  SFREE(nelemloadref);if(*nam>0){SFREE(iamloadref);};
1523  SFREE(sideloadref);
1524  }
1525 
1526  /* setting the velocity to zero at the end of a quasistatic or stationary
1527  step */
1528 
1529  if(abs(*nmethod)==1){
1530  for(k=0;k<mt**nk;++k){veold[k]=0.;}
1531  }
1532 
1533  /* updating the loading at the end of the step;
1534  important in case the amplitude at the end of the step
1535  is not equal to one */
1536 
1537  for(k=0;k<*nboun;++k){
1538 
1539  /* thermal boundary conditions are updated only if the
1540  step was thermal or thermomechanical */
1541 
1542  if(ndirboun[k]==0){
1543  if(*ithermal<2) continue;
1544 
1545  /* mechanical boundary conditions are updated only
1546  if the step was not thermal or the node is a
1547  network node */
1548 
1549  }else if((ndirboun[k]>0)&&(ndirboun[k]<4)){
1550  node=nodeboun[k];
1551  FORTRAN(nident,(itg,&node,&ntg,&id));
1552  networknode=0;
1553  if(id>0){
1554  if(itg[id-1]==node) networknode=1;
1555  }
1556  if((*ithermal==2)&&(networknode==0)) continue;
1557  }
1558  xbounold[k]=xbounact[k];
1559  }
1560  for(k=0;k<*nforc;++k){xforcold[k]=xforcact[k];}
1561  for(k=0;k<2**nload;++k){xloadold[k]=xloadact[k];}
1562  for(k=0;k<7**nbody;k=k+7){xbodyold[k]=xbodyact[k];}
1563  if(*ithermal==1){
1564  for(k=0;k<*nk;++k){t1old[k]=t1act[k];}
1565  for(k=0;k<*nk;++k){vold[mt*k]=t1act[k];}
1566  }
1567  else if(*ithermal>1){
1568  for(k=0;k<*nk;++k){t1[k]=vold[mt*k];}
1569  if(*ithermal>=3){
1570  for(k=0;k<*nk;++k){t1old[k]=t1act[k];}
1571  }
1572  }
1573 
1574  qaold[0]=qa[0];
1575  qaold[1]=qa[1];
1576 
1577  SFREE(f);
1578  SFREE(b);
1579  SFREE(xbounact);SFREE(xforcact);SFREE(xloadact);SFREE(xbodyact);
1580  if(*nbody>0) SFREE(ipobody);
1581  SFREE(fext);SFREE(ampli);SFREE(xbounini);SFREE(xstiff);
1582  if((*ithermal==1)||(*ithermal>=3)){SFREE(t1act);SFREE(t1ini);}
1583 
1584  if(*ithermal>1){
1585  SFREE(itg);SFREE(ieg);SFREE(kontri);SFREE(nloadtr);
1586  SFREE(nactdog);SFREE(nacteq);SFREE(ineighe);
1587  SFREE(tarea);SFREE(tenv);SFREE(fenv);SFREE(qfx);
1588  SFREE(erad);SFREE(ac);SFREE(bc);SFREE(ipiv);
1589  SFREE(bcr);SFREE(ipivr);SFREE(adview);SFREE(auview);SFREE(adrad);
1590  SFREE(aurad);SFREE(irowrad);SFREE(jqrad);SFREE(icolrad);
1591  if((*mcs>0)&&(ntr>0)){SFREE(inocs);}
1592  if((*network>0)||(ntg>0)){SFREE(iponoel);SFREE(inoel);}
1593  }
1594 
1595  SFREE(fini);
1596  if(*nmethod==4){
1597  SFREE(aux2);SFREE(fextini);SFREE(veini);
1598  SFREE(adb);SFREE(aub);SFREE(cv);SFREE(cvini);
1599  }
1600 
1601  SFREE(aux);SFREE(iaux);SFREE(vini);SFREE(h0ref);SFREE(h0);SFREE(inomat);
1602 
1603  /* reset icascade */
1604 
1605  if(icascade==1){icascade=0;}
1606 
1607  mpcinfo[0]=memmpc_;mpcinfo[1]=mpcfree;mpcinfo[2]=icascade;
1608  mpcinfo[3]=maxlenmpc;
1609 
1610  *icolp=icol;*irowp=irow;*cop=co;*voldp=vold;
1611 
1612  *ipompcp=ipompc;*labmpcp=labmpc;*ikmpcp=ikmpc;*ilmpcp=ilmpc;
1613  *fmpcp=fmpc;*nodempcp=nodempc;*coefmpcp=coefmpc;
1614 
1615  *ipkonp=ipkon;*lakonp=lakon;*konp=kon;*ielorienp=ielorien;
1616  *ielmatp=ielmat;*enerp=ener;*xstatep=xstate;
1617 
1618  *setp=set;*istartsetp=istartset;*iendsetp=iendset;*ialsetp=ialset;
1619  *tiesetp=tieset;*tietolp=tietol;
1620 
1621  *nelemloadp=nelemload;*iamloadp=iamload;
1622  *sideloadp=sideload;
1623 
1624  (*tmin)*=(*tper);
1625  (*tmax)*=(*tper);
1626 
1627  SFREE(nactdofinv);
1628 
1629  if(*nmethod==1){
1630  *nmethod=8;
1631  }else if(*nmethod==4){
1632  *nmethod=9;
1633  }else if(*nmethod==2){
1634  *nmethod=10;
1635  }
1636 
1637  (*ttime)+=(*tper);
1638 
1639  return;
1640 }
subroutine checktime(itpamp, namta, tinc, ttime, amta, tmin, inext, itp, istep, tper)
Definition: checktime.f:21
#define ITGFORMAT
Definition: CalculiX.h:52
void pardiso_main(double *ad, double *au, double *adb, double *aub, double *sigma, double *b, ITG *icol, ITG *irow, ITG *neq, ITG *nzs, ITG *symmetryflag, ITG *inputformat, ITG *jq, ITG *nzs3)
subroutine assigndomtonodes(ne, lakon, ipkon, kon, ielmat, inomat, elcon, ncmat_, ntmat_, mi, ne2)
Definition: assigndomtonodes.f:21
subroutine createtiedsurfs(nodface, ipoface, set, istartset, iendset, ialset, tieset, inomat, ne, ipkon, lakon, kon, ntie, tietol, nalset, nk, nset, iactive)
Definition: createtiedsurfs.f:22
void frd(double *co, ITG *nk, ITG *kon, ITG *ipkon, char *lakon, ITG *ne0, double *v, double *stn, ITG *inum, ITG *nmethod, ITG *kode, char *filab, double *een, double *t1, double *fn, double *time, double *epn, ITG *ielmat, char *matname, double *enern, double *xstaten, ITG *nstate_, ITG *istep, ITG *iinc, ITG *ithermal, double *qfn, ITG *mode, ITG *noddiam, double *trab, ITG *inotr, ITG *ntrans, double *orab, ITG *ielorien, ITG *norien, char *description, ITG *ipneigh, ITG *neigh, ITG *mi, double *stx, double *vr, double *vi, double *stnr, double *stni, double *vmax, double *stnmax, ITG *ngraph, double *veold, double *ener, ITG *ne, double *cs, char *set, ITG *nset, ITG *istartset, ITG *iendset, ITG *ialset, double *eenmax, double *fnr, double *fni, double *emn, double *thicke, char *jobnamec, char *output, double *qfx, double *cdn, ITG *mortar, double *cdnr, double *cdni, ITG *nmat)
Definition: frd.c:32
void mafillsmmain(double *co, ITG *nk, ITG *kon, ITG *ipkon, char *lakon, ITG *ne, ITG *nodeboun, ITG *ndirboun, double *xboun, ITG *nboun, ITG *ipompc, ITG *nodempc, double *coefmpc, ITG *nmpc, ITG *nodeforc, ITG *ndirforc, double *xforc, ITG *nforc, ITG *nelemload, char *sideload, double *xload, ITG *nload, double *xbody, ITG *ipobody, ITG *nbody, double *cgr, double *ad, double *au, double *bb, ITG *nactdof, ITG *icol, ITG *jq, ITG *irow, ITG *neq, ITG *nzl, ITG *nmethod, ITG *ikmpc, ITG *ilmpc, ITG *ikboun, ITG *ilboun, double *elcon, ITG *nelcon, double *rhcon, ITG *nrhcon, double *alcon, ITG *nalcon, double *alzero, ITG *ielmat, ITG *ielorien, ITG *norien, double *orab, ITG *ntmat_, double *t0, double *t1, ITG *ithermal, double *prestr, ITG *iprestr, double *vold, ITG *iperturb, double *sti, ITG *nzs, double *stx, double *adb, double *aub, ITG *iexpl, double *plicon, ITG *nplicon, double *plkcon, ITG *nplkcon, double *xstiff, ITG *npmat_, double *dtime, char *matname, ITG *mi, ITG *ncmat_, ITG *mass, ITG *stiffness, ITG *buckling, ITG *rhs, ITG *intscheme, double *physcon, double *shcon, ITG *nshcon, double *cocon, ITG *ncocon, double *ttime, double *time, ITG *istep, ITG *kinc, ITG *coriolis, ITG *ibody, double *xloadold, double *reltime, double *veold, double *springarea, ITG *nstate_, double *xstateini, double *xstate, double *thicke, ITG *integerglob, double *doubleglob, char *tieset, ITG *istartset, ITG *iendset, ITG *ialset, ITG *ntie, ITG *nasym, double *pslavsurf, double *pmastsurf, ITG *mortar, double *clearini, ITG *ielprop, double *prop, ITG *ne0, double *fnext, ITG *kscale, ITG *iponoel, ITG *inoel, ITG *network)
Definition: mafillsmmain.c:47
void preiter(double *ad, double **aup, double *b, ITG **icolp, ITG **irowp, ITG *neq, ITG *nzs, ITG *isolver, ITG *iperturb)
Definition: preiter.c:23
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))
void prediction_em(double *uam, ITG *nmethod, double *bet, double *gam, double *dtime, ITG *ithermal, ITG *nk, double *veold, double *v, ITG *iinc, ITG *idiscon, double *vold, ITG *nactdof, ITG *mi)
Definition: prediction_em.c:33
void calcresidual_em(ITG *nmethod, ITG *neq, double *b, double *fext, double *f, ITG *iexpl, ITG *nactdof, double *aux1, double *aux2, double *vold, double *vini, double *dtime, double *accold, ITG *nk, double *adb, double *aub, ITG *jq, ITG *irow, ITG *nzl, double *alpha, double *fextini, double *fini, ITG *islavnode, ITG *nslavnode, ITG *mortar, ITG *ntie, double *f_cm, double *f_cs, ITG *mi, ITG *nzs, ITG *nasym, ITG *ithermal)
Definition: calcresidual_em.c:33
ITG strcmp1(const char *s1, const char *s2)
Definition: strcmp1.c:24
void tiedcontact(ITG *ntie, char *tieset, ITG *nset, char *set, ITG *istartset, ITG *iendset, ITG *ialset, char *lakon, ITG *ipkon, ITG *kon, double *tietol, ITG *nmpc, ITG *mpcfree, ITG *memmpc_, ITG **ipompcp, char **labmpcp, ITG **ikmpcp, ITG **ilmpcp, double **fmpcp, ITG **nodempcp, double **coefmpcp, ITG *ithermal, double *co, double *vold, ITG *cfd, ITG *nmpc_, ITG *mi, ITG *nk, ITG *istep, ITG *ikboun, ITG *nboun, char *kind1, char *kind2)
Definition: tiedcontact.c:23
void mastructrad(ITG *ntr, ITG *nloadtr, char *sideload, ITG *ipointerrad, ITG **mast1radp, ITG **irowradp, ITG *nzsrad, ITG *jqrad, ITG *icolrad)
Definition: mastructrad.c:24
subroutine preconditioning(ad, au, b, neq, irow, jq, adaux)
Definition: preconditioning.f:22
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
void spooles(double *ad, double *au, double *adb, double *aub, double *sigma, double *b, ITG *icol, ITG *irow, ITG *neq, ITG *nzs, ITG *symmtryflag, ITG *inputformat, ITG *nzs3)
subroutine stop()
Definition: stop.f:20
void radflowload(ITG *itg, ITG *ieg, ITG *ntg, ITG *ntr, double *adrad, double *aurad, double *bcr, ITG *ipivr, double *ac, double *bc, ITG *nload, char *sideload, ITG *nelemload, double *xloadact, char *lakon, ITG *ipiv, ITG *ntmat_, double *vold, double *shcon, ITG *nshcon, ITG *ipkon, ITG *kon, double *co, ITG *kontri, ITG *ntri, ITG *nloadtr, double *tarea, double *tenv, double *physcon, double *erad, double **adviewp, double **auviewp, ITG *nflow, ITG *ikboun, double *xboun, ITG *nboun, ITG *ithermal, ITG *iinc, ITG *iit, double *cs, ITG *mcs, ITG *inocs, ITG *ntrit, ITG *nk, double *fenv, ITG *istep, double *dtime, double *ttime, double *time, ITG *ilboun, ITG *ikforc, ITG *ilforc, double *xforc, ITG *nforc, double *cam, ITG *ielmat, ITG *nteq, double *prop, ITG *ielprop, ITG *nactdog, ITG *nacteq, ITG *nodeboun, ITG *ndirboun, ITG *network, double *rhcon, ITG *nrhcon, ITG *ipobody, ITG *ibody, double *xbody, ITG *nbody, ITG *iviewfile, char *jobnamef, double *ctrl, double *xloadold, double *reltime, ITG *nmethod, char *set, ITG *mi, ITG *istartset, ITG *iendset, ITG *ialset, ITG *nset, ITG *ineighe, ITG *nmpc, ITG *nodempc, ITG *ipompc, double *coefmpc, char *labmpc, ITG *iemchange, ITG *nam, ITG *iamload, ITG *jqrad, ITG *irowrad, ITG *nzsrad, ITG *icolrad, ITG *ne, ITG *iaxial, double *qa, double *cocon, ITG *ncocon, ITG *iponoel, ITG *inoel, ITG *nprop, char *amname, ITG *namta, double *amta)
Definition: radflowload.c:45
void tau(double *ad, double **aup, double *adb, double *aubp, double *sigma, double *b, ITG *icol, ITG **irowp, ITG *neq, ITG *nzs)
subroutine createinterfacempcs(imastnode, xmastnor, nmastnode, ikmpc, ilmpc, nmpc, ipompc, nodempc, coefmpc, labmpc, mpcfree, ikboun, nboun)
Definition: createinterfacempcs.f:22
void radcyc(ITG *nk, ITG *kon, ITG *ipkon, char *lakon, ITG *ne, double *cs, ITG *mcs, ITG *nkon, ITG *ialset, ITG *istartset, ITG *iendset, ITG **kontrip, ITG *ntri, double **cop, double **voldp, ITG *ntrit, ITG *inocs, ITG *mi)
Definition: radcyc.c:24
void sgi_main(double *ad, double *au, double *adb, double *aub, double *sigma, double *b, ITG *icol, ITG *irow, ITG *neq, ITG *nzs, ITG token)
#define RENEW(a, b, c)
Definition: CalculiX.h:40
subroutine tempload_em(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, h0scale, inomat, ipobody, iponoel, inoel)
Definition: tempload_em.f:29
#define SFREE(a)
Definition: CalculiX.h:41
subroutine gennactdofinv(nactdof, nactdofinv, nk, mi, nodorig, ipkon, lakon, kon, ne)
Definition: gennactdofinv.f:21
subroutine calch0interface(nmpc, ipompc, nodempc, coefmpc, h0)
Definition: calch0interface.f:27
void frdcyc(double *co, ITG *nk, ITG *kon, ITG *ipkon, char *lakon, ITG *ne, double *v, double *stn, ITG *inum, ITG *nmethod, ITG *kode, char *filab, double *een, double *t1, double *fn, double *time, double *epn, ITG *ielmat, char *matname, double *cs, ITG *mcs, ITG *nkon, double *enern, double *xstaten, ITG *nstate_, ITG *istep, ITG *iinc, ITG *iperturb, double *ener, ITG *mi, char *output, ITG *ithermal, double *qfn, ITG *ialset, ITG *istartset, ITG *iendset, double *trab, ITG *inotr, ITG *ntrans, double *orab, ITG *ielorien, ITG *norien, double *sti, double *veold, ITG *noddiam, char *set, ITG *nset, double *emn, double *thicke, char *jobnamec, ITG *ne0, double *cdn, ITG *mortar, ITG *nmat, double *qfx)
Definition: frdcyc.c:24
void resultsinduction(double *co, ITG *nk, ITG *kon, ITG *ipkon, char *lakon, ITG *ne, double *v, double *stn, ITG *inum, double *elcon, ITG *nelcon, double *rhcon, ITG *nrhcon, double *alcon, ITG *nalcon, double *alzero, ITG *ielmat, ITG *ielorien, ITG *norien, double *orab, ITG *ntmat_, double *t0, double *t1, ITG *ithermal, double *prestr, ITG *iprestr, char *filab, double *eme, double *emn, double *een, ITG *iperturb, double *f, double *fn, ITG *nactdof, ITG *iout, double *qa, double *vold, double *b, ITG *nodeboun, ITG *ndirboun, double *xboun, ITG *nboun, ITG *ipompc, ITG *nodempc, double *coefmpc, char *labmpc, ITG *nmpc, ITG *nmethod, double *vmax, ITG *neq, double *veold, double *accold, double *beta, double *gamma, double *dtime, double *time, double *ttime, double *plicon, ITG *nplicon, double *plkcon, ITG *nplkcon, double *xstateini, double *xstiff, double *xstate, ITG *npmat_, double *epl, char *matname, ITG *mi, ITG *ielas, ITG *icmd, ITG *ncmat_, ITG *nstate_, double *sti, double *vini, ITG *ikboun, ITG *ilboun, double *ener, double *enern, double *emeini, double *xstaten, double *eei, double *enerini, double *cocon, ITG *ncocon, char *set, ITG *nset, ITG *istartset, ITG *iendset, ITG *ialset, ITG *nprint, char *prlab, char *prset, double *qfx, double *qfn, double *trab, ITG *inotr, ITG *ntrans, double *fmpc, ITG *nelemload, ITG *nload, ITG *ikmpc, ITG *ilmpc, ITG *istep, ITG *iinc, double *springarea, double *reltime, ITG *ne0, double *xforc, ITG *nforc, double *thicke, double *shcon, ITG *nshcon, char *sideload, double *xload, double *xloadold, ITG *icfd, ITG *inomat, double *h0, ITG *islavnode, ITG *nslavnode, ITG *ntie, ITG *ielprop, double *prop, ITG *iactive, double *energyini, double *energy, ITG *iponoel, ITG *inoel, char *orname, ITG *network, ITG *ipobody, double *xbody, ITG *ibody)
Definition: resultsinduction.c:42
void remastruct(ITG *ipompc, double **coefmpcp, ITG **nodempcp, ITG *nmpc, ITG *mpcfree, ITG *nodeboun, ITG *ndirboun, ITG *nboun, ITG *ikmpc, ITG *ilmpc, ITG *ikboun, ITG *ilboun, char *labmpc, ITG *nk, ITG *memmpc_, ITG *icascade, ITG *maxlenmpc, ITG *kon, ITG *ipkon, char *lakon, ITG *ne, ITG *nactdof, ITG *icol, ITG *jq, ITG **irowp, ITG *isolver, ITG *neq, ITG *nzs, ITG *nmethod, double **fp, double **fextp, double **bp, double **aux2p, double **finip, double **fextinip, double **adbp, double **aubp, ITG *ithermal, ITG *iperturb, ITG *mass, ITG *mi, ITG *iexpl, ITG *mortar, char *typeboun, double **cvp, double **cvinip, ITG *iit, ITG *network)
Definition: remastruct.c:24
static double * adview
Definition: radflowload.c:42
subroutine nident(x, px, n, id)
Definition: nident.f:26
real *8 function f_cm(x, phi, lambda1, zk0, Pup, Tup, rurd, xflow, kup)
Definition: moehring.f:582
subroutine envtemp(itg, ieg, ntg, ntr, sideload, nelemload, ipkon, kon, lakon, ielmat, ne, nload, kontri, ntri, nloadtr, nflow, ndirboun, nactdog, nodeboun, nacteq, nboun, ielprop, prop, nteq, v, network, physcon, shcon, ntmat_, co, vold, set, nshcon, rhcon, nrhcon, mi, nmpc, nodempc, ipompc, labmpc, ikboun, nasym, ttime, time, iaxial)
Definition: envtemp.f:25
void remastructem(ITG *ipompc, double **coefmpcp, ITG **nodempcp, ITG *nmpc, ITG *mpcfree, ITG *nodeboun, ITG *ndirboun, ITG *nboun, ITG *ikmpc, ITG *ilmpc, ITG *ikboun, ITG *ilboun, char *labmpc, ITG *nk, ITG *memmpc_, ITG *icascade, ITG *maxlenmpc, ITG *kon, ITG *ipkon, char *lakon, ITG *ne, ITG *nactdof, ITG *icol, ITG *jq, ITG **irowp, ITG *isolver, ITG *neq, ITG *nzs, ITG *nmethod, double **fp, double **fextp, double **bp, double **aux2p, double **finip, double **fextinip, double **adbp, double **aubp, ITG *ithermal, ITG *iperturb, ITG *mass, ITG *mi, ITG *ielmat, double *elcon, ITG *ncmat_, ITG *ntmat_, ITG *inomat, ITG *network)
Definition: remastructem.c:24
subroutine createinum(ipkon, inum, kon, lakon, nk, ne, cflag, nelemload, nload, nodeboun, nboun, ndirboun, ithermal, co, vold, mi, ielmat)
Definition: createinum.f:21
subroutine jouleheating(ipkon, lakon, kon, co, elcon, nelcon, mi, ne, sti, ielmat, nelemload, sideload, xload, nload, nload_, iamload, nam, idefload, ncmat_, ntmat_, alcon, nalcon, ithermal, vold, t1)
Definition: jouleheating.f:23
void biosav(ITG *ipkon, ITG *kon, char *lakon, ITG *ne, double *co, double *qfx, double *h0, ITG *mi, ITG *inomat, ITG *nk)
Definition: biosav.c:31
subroutine networkelementpernode(iponoel, inoel, lakon, ipkon, kon, inoelsize, nflow, ieg, ne, network)
Definition: networkelementpernode.f:21
subroutine generateeminterfaces(istartset, iendset, ialset, iactive, ipkon, lakon, kon, ikmpc, nmpc, nafaces)
Definition: generateeminterfaces.f:21
void checkconvergence(double *co, ITG *nk, ITG *kon, ITG *ipkon, char *lakon, ITG *ne, double *stn, ITG *nmethod, ITG *kode, char *filab, double *een, double *t1act, double *time, double *epn, ITG *ielmat, char *matname, double *enern, double *xstaten, ITG *nstate_, ITG *istep, ITG *iinc, ITG *iperturb, double *ener, ITG *mi, char *output, ITG *ithermal, double *qfn, ITG *mode, ITG *noddiam, double *trab, ITG *inotr, ITG *ntrans, double *orab, ITG *ielorien, ITG *norien, char *description, double *sti, ITG *icutb, ITG *iit, double *dtime, double *qa, double *vold, double *qam, double *ram1, double *ram2, double *ram, double *cam, double *uam, ITG *ntg, double *ttime, ITG *icntrl, double *theta, double *dtheta, double *veold, double *vini, ITG *idrct, double *tper, ITG *istab, double *tmax, ITG *nactdof, double *b, double *tmin, double *ctrl, double *amta, ITG *namta, ITG *itpamp, ITG *inext, double *dthetaref, ITG *itp, ITG *jprint, ITG *jout, ITG *uncoupled, double *t1, ITG *iitterm, ITG *nelemload, ITG *nload, ITG *nodeboun, ITG *nboun, ITG *itg, ITG *ndirboun, double *deltmx, ITG *iflagact, char *set, ITG *nset, ITG *istartset, ITG *iendset, ITG *ialset, double *emn, double *thicke, char *jobnamec, ITG *mortar, ITG *nmat, ITG *ielprop, double *prop, ITG *ialeatoric, ITG *kscale, double *energy, double *allwk, double *energyref, double *emax, double *enres, double *enetoll, double *energyini, double *allwkini, double *temax, double *reswk, ITG *ne0, ITG *neini, double *dampwk, double *dampwkini, double *energystartstep)
Definition: checkconvergence.c:34
#define ITG
Definition: CalculiX.h:51
void results(double *co, ITG *nk, ITG *kon, ITG *ipkon, char *lakon, ITG *ne, double *v, double *stn, ITG *inum, double *stx, double *elcon, ITG *nelcon, double *rhcon, ITG *nrhcon, double *alcon, ITG *nalcon, double *alzero, ITG *ielmat, ITG *ielorien, ITG *norien, double *orab, ITG *ntmat_, double *t0, double *t1, ITG *ithermal, double *prestr, ITG *iprestr, char *filab, double *eme, double *emn, double *een, ITG *iperturb, double *f, double *fn, ITG *nactdof, ITG *iout, double *qa, double *vold, double *b, ITG *nodeboun, ITG *ndirboun, double *xboun, ITG *nboun, ITG *ipompc, ITG *nodempc, double *coefmpc, char *labmpc, ITG *nmpc, ITG *nmethod, double *vmax, ITG *neq, double *veold, double *accold, double *beta, double *gamma, double *dtime, double *time, double *ttime, double *plicon, ITG *nplicon, double *plkcon, ITG *nplkcon, double *xstateini, double *xstiff, double *xstate, ITG *npmat_, double *epl, char *matname, ITG *mi, ITG *ielas, ITG *icmd, ITG *ncmat_, ITG *nstate_, double *stiini, double *vini, ITG *ikboun, ITG *ilboun, double *ener, double *enern, double *emeini, double *xstaten, double *eei, double *enerini, double *cocon, ITG *ncocon, char *set, ITG *nset, ITG *istartset, ITG *iendset, ITG *ialset, ITG *nprint, char *prlab, char *prset, double *qfx, double *qfn, double *trab, ITG *inotr, ITG *ntrans, double *fmpc, ITG *nelemload, ITG *nload, ITG *ikmpc, ITG *ilmpc, ITG *istep, ITG *iinc, double *springarea, double *reltime, ITG *ne0, double *xforc, ITG *nforc, double *thicke, double *shcon, ITG *nshcon, char *sideload, double *xload, double *xloadold, ITG *icfd, ITG *inomat, double *pslavsurf, double *pmastsurf, ITG *mortar, ITG *islavact, double *cdn, ITG *islavnode, ITG *nslavnode, ITG *ntie, double *clearini, ITG *islavsurf, ITG *ielprop, double *prop, double *energyini, double *energy, ITG *kscale, ITG *iponoel, ITG *inoel, ITG *nener, char *orname, ITG *network, ITG *ipobody, double *xbodyact, ITG *ibody)
Definition: results.c:42
subroutine mafillem(co, nk, kon, ipkon, lakon, ne, nodeboun, ndirboun, xboun, nboun, ipompc, nodempc, coefmpc, nmpc, nodeforc, ndirforc, xforc, nforc, nelemload, sideload, xload, nload, xbody, ipobody, nbody, cgr, ad, au, fext, nactdof, icol, jq, irow, neq, nzl, nmethod, ikmpc, ilmpc, ikboun, ilboun, elcon, nelcon, rhcon, nrhcon, alcon, nalcon, alzero, ielmat, ielorien, norien, orab, ntmat_, t0, t1, ithermal, prestr, iprestr, vold, iperturb, sti, nzs, stx, adb, aub, iexpl, plicon, nplicon, plkcon, nplkcon, xstiff, npmat_, dtime, matname, mi, ncmat_, mass, stiffness, buckling, rhsi, intscheme, physcon, shcon, nshcon, cocon, ncocon, ttime, time, istep, iinc, coriolis, ibody, xloadold, reltime, veold, springarea, nstate_, xstateini, xstate, thicke, integerglob, doubleglob, tieset, istartset, iendset, ialset, ntie, nasym, iactive, h0, pslavsurf, pmastsurf, mortar, clearini, ielprop, prop, iponoel, inoel, network)
Definition: mafillem.f:36
#define NNEW(a, b, c)
Definition: CalculiX.h:39
subroutine bodyforce(cbody, ibody, ipobody, nbody, set, istartset, iendset, ialset, inewton, nset, ifreebody, k)
Definition: bodyforce.f:21
static double * auview
Definition: radflowload.c:42
subroutine normalsoninterface(istartset, iendset, ialset, imast, ipkon, kon, lakon, imastnode, nmastnode, xmastnor, co)
Definition: normalsoninterface.f:22
Hosted by OpenAircraft.com, (Michigan UAV, LLC)