CalculiX  2.13
A Free Software Three-Dimensional Structural Finite Element Program
resultsinduction.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 resultsinduction.c:

Go to the source code of this file.

Functions

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 *cam, ITG *neq, double *veold, double *accold, double *bet, double *gam, double *dtime, double *time, double *ttime, double *plicon, ITG *nplicon, double *plkcon, ITG *nplkcon, double *xstateini, double *xstiff, double *xstate, ITG *npmat_, double *epn, 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)
 
void * resultsemmt (ITG *i)
 
void * resultsthermemmt (ITG *i)
 

Variables

static char * lakon1
 
static char * matname1
 
static char * sideload1
 
static ITGkon1
 
static ITGipkon1
 
static ITGne1
 
static ITGnelcon1
 
static ITGnrhcon1
 
static ITGnalcon1
 
static ITGielmat1
 
static ITGielorien1
 
static ITGnorien1
 
static ITGntmat1_
 
static ITGithermal1
 
static ITGiperturb1
 
static ITGiout1
 
static ITGnmethod1
 
static ITGnplkcon1
 
static ITGnpmat1_
 
static ITGmi1
 
static ITGncmat1_
 
static ITGnstate1_
 
static ITGielprop1
 
static ITGinoel1
 
static ITGistep1
 
static ITGiinc1
 
static ITG calcul_fn1
 
static ITG calcul_qa1
 
static ITGnplicon1
 
static ITGiponoel1
 
static ITGnal =NULL
 
static ITGipompc1
 
static ITGnodempc1
 
static ITGnmpc1
 
static ITGncocon1
 
static ITGikmpc1
 
static ITGilmpc1
 
static ITG num_cpus
 
static ITG mt1
 
static ITGnk1
 
static ITGnshcon1
 
static ITGnelemload1
 
static ITGnload1
 
static ITG mortar1
 
static ITGistartset1
 
static ITGiendset1
 
static ITGialset1
 
static ITGiactive1
 
static ITGnetwork1
 
static ITGipobody1
 
static ITGibody1
 
static double * co1
 
static double * v1
 
static double * elcon1
 
static double * rhcon1
 
static double * alcon1
 
static double * orab1
 
static double * t01
 
static double * fn1 =NULL
 
static double * qa1 =NULL
 
static double * vold1
 
static double * dtime1
 
static double * time1
 
static double * prop1
 
static double * ttime1
 
static double * plkcon1
 
static double * xstateini1
 
static double * xstiff1
 
static double * xstate1
 
static double * sti1
 
static double * springarea1
 
static double * reltime1
 
static double * coefmpc1
 
static double * vini1
 
static double * cocon1
 
static double * qfx1
 
static double * shcon1
 
static double * xload1
 
static double * plicon1
 
static double * xloadold1
 
static double * h01
 
static double * pslavsurf1
 
static double * pmastsurf1
 
static double * clearini1
 
static double * xbody1
 

Function Documentation

◆ resultsemmt()

void* resultsemmt ( ITG i)
311  {
312 
313  ITG nea,neb,nedelta;
314 
315  nedelta=(ITG)floor(*ne1/(double)num_cpus);
316  nea=*i*nedelta+1;
317  neb=(*i+1)*nedelta;
318 // next line! -> all parallel sections
319  if((*i==num_cpus-1)&&(neb<*ne1)) neb=*ne1;
320 
324 
325  return NULL;
326 }
static ITG * mi1
Definition: resultsinduction.c:27
static double * v1
Definition: resultsinduction.c:35
static double * sti1
Definition: resultsinduction.c:36
static ITG num_cpus
Definition: resultsinduction.c:31
static ITG * istartset1
Definition: resultsinduction.c:31
static double * elcon1
Definition: resultsinduction.c:35
static double * vini1
Definition: resultsinduction.c:36
static char * matname1
Definition: resultsinduction.c:25
static ITG * nelcon1
Definition: resultsinduction.c:27
void FORTRAN(actideacti,(char *set, ITG *nset, ITG *istartset, ITG *iendset, ITG *ialset, char *objectset, ITG *ipkon, ITG *ibject, ITG *ne))
static ITG * kon1
Definition: resultsinduction.c:27
static double * alcon1
Definition: resultsinduction.c:35
static ITG * ipkon1
Definition: resultsinduction.c:27
static ITG * nalcon1
Definition: resultsinduction.c:27
subroutine resultsem(co, kon, ipkon, lakon, v, elcon, nelcon, ielmat, ntmat_, vini, dtime, matname, mi, ncmat_, nea, neb, sti, alcon, nalcon, h0, istartset, iendset, ialset, iactive, fn)
Definition: resultsem.f:22
static double * fn1
Definition: resultsinduction.c:36
static ITG * ielmat1
Definition: resultsinduction.c:27
static ITG * iendset1
Definition: resultsinduction.c:31
static char * lakon1
Definition: resultsinduction.c:25
#define ITG
Definition: CalculiX.h:51
static ITG * ialset1
Definition: resultsinduction.c:31
static ITG * ncmat1_
Definition: resultsinduction.c:27
static double * co1
Definition: resultsinduction.c:35
static ITG * ne1
Definition: resultsinduction.c:27
static double * dtime1
Definition: resultsinduction.c:36
static ITG * ntmat1_
Definition: resultsinduction.c:27
static ITG * iactive1
Definition: resultsinduction.c:31
static double * h01
Definition: resultsinduction.c:36

◆ resultsinduction()

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 *  cam,
ITG neq,
double *  veold,
double *  accold,
double *  bet,
double *  gam,
double *  dtime,
double *  time,
double *  ttime,
double *  plicon,
ITG nplicon,
double *  plkcon,
ITG nplkcon,
double *  xstateini,
double *  xstiff,
double *  xstate,
ITG npmat_,
double *  epn,
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 
)
73  {
74 
75  /* variables for multithreading procedure */
76 
77  char *env,*envloc,*envsys;
78 
79  ITG intpointvarm,calcul_fn,calcul_f,calcul_qa,calcul_cauchy,nener,ikin,
80  intpointvart,mt=mi[1]+1,i,j,*ithread=NULL,*islavsurf=NULL,
81  sys_cpus,mortar=0,*islavact=NULL;
82 
83  double *pmastsurf=NULL,*clearini=NULL,*pslavsurf=NULL,*cdn=NULL;
84 
85  /*
86 
87  calculating integration point values (strains, stresses,
88  heat fluxes, material tangent matrices and nodal forces)
89 
90  storing the nodal and integration point results in the
91  .dat file
92 
93  iout=-2: v is assumed to be known and is used to
94  calculate strains, stresses..., no result output
95  corresponds to iout=-1 with in addition the
96  calculation of the internal energy density
97  iout=-1: v is assumed to be known and is used to
98  calculate strains, stresses..., no result output;
99  is used to take changes in SPC's and MPC's at the
100  start of a new increment or iteration into account
101  iout=0: v is calculated from the system solution
102  and strains, stresses.. are calculated, no result output
103  iout=1: v is calculated from the system solution and strains,
104  stresses.. are calculated, requested results output
105  iout=2: v is assumed to be known and is used to
106  calculate strains, stresses..., requested results output */
107 
108  num_cpus=0;
109  sys_cpus=0;
110 
111  /* explicit user declaration prevails */
112 
113  envsys=getenv("NUMBER_OF_CPUS");
114  if(envsys){
115  sys_cpus=atoi(envsys);
116  if(sys_cpus<0) sys_cpus=0;
117  }
118 
119  /* automatic detection of available number of processors */
120 
121  if(sys_cpus==0){
122  sys_cpus = getSystemCPUs();
123  if(sys_cpus<1) sys_cpus=1;
124  }
125 
126  /* local declaration prevails, if strictly positive */
127 
128  envloc = getenv("CCX_NPROC_RESULTS");
129  if(envloc){
130  num_cpus=atoi(envloc);
131  if(num_cpus<0){
132  num_cpus=0;
133  }else if(num_cpus>sys_cpus){
134  num_cpus=sys_cpus;
135  }
136 
137  }
138 
139  /* else global declaration, if any, applies */
140 
141  env = getenv("OMP_NUM_THREADS");
142  if(num_cpus==0){
143  if (env)
144  num_cpus = atoi(env);
145  if (num_cpus < 1) {
146  num_cpus=1;
147  }else if(num_cpus>sys_cpus){
148  num_cpus=sys_cpus;
149  }
150  }
151 
152 // next line is to be inserted in a similar way for all other paralell parts
153 
154  if(*ne<num_cpus) num_cpus=*ne;
155 
156  pthread_t tid[num_cpus];
157 
158  /* 1. nodewise storage of the primary variables
159  2. determination which derived variables have to be calculated */
160 
161  FORTRAN(resultsini_em,(nk,v,ithermal,filab,iperturb,f,fn,
162  nactdof,iout,qa,b,nodeboun,ndirboun,
163  xboun,nboun,ipompc,nodempc,coefmpc,labmpc,nmpc,nmethod,cam,neq,
164  veold,dtime,mi,vini,nprint,prlab,
165  &intpointvarm,&calcul_fn,&calcul_f,&calcul_qa,&calcul_cauchy,&nener,
166  &ikin,&intpointvart,xforc,nforc));
167 
168  /* electromagnetic calculation is linear: should not be taken
169  into account in the convergence check (only thermal part
170  is taken into account) */
171 
172  cam[0]=0.;
173 
174  /* next statement allows for storing the displacements in each
175  iteration: for debugging purposes */
176 
177  if((strcmp1(&filab[3],"I")==0)&&(*iout==0)){
178  FORTRAN(frditeration,(co,nk,kon,ipkon,lakon,ne,v,
179  ttime,ielmat,matname,mi,istep,iinc,ithermal));
180  }
181 
182  /* calculating the stresses and material tangent at the
183  integration points; calculating the internal forces */
184 
185  if(((ithermal[0]<=1)||(ithermal[0]>=3))&&(intpointvarm==1)){
186 
187  co1=co;kon1=kon;ipkon1=ipkon;lakon1=lakon;v1=v;elcon1=elcon;
188  nelcon1=nelcon;ielmat1=ielmat;ntmat1_=ntmat_;vini1=vini;dtime1=dtime;
189  matname1=matname;mi1=mi;ncmat1_=ncmat_;sti1=sti;alcon1=alcon;
190  nalcon1=nalcon;h01=h0;ne1=ne;istartset1=istartset;iendset1=iendset;
191  ialset1=ialset;iactive1=iactive;fn1=fn;
192 
193  /* calculating the magnetic field */
194 
195  if(((*nmethod!=4)&&(*nmethod!=5))||(iperturb[0]>1)){
196  printf(" Using up to %" ITGFORMAT " cpu(s) for the magnetic field calculation.\n\n", num_cpus);
197  }
198 
199  /* create threads and wait */
200 
201  NNEW(ithread,ITG,num_cpus);
202  for(i=0; i<num_cpus; i++) {
203  ithread[i]=i;
204  pthread_create(&tid[i], NULL, (void *)resultsemmt, (void *)&ithread[i]);
205  }
206  for(i=0; i<num_cpus; i++) pthread_join(tid[i], NULL);
207  SFREE(ithread);
208 
209  qa[0]=0.;
210  }
211 
212  /* calculating the thermal flux and material tangent at the
213  integration points; calculating the internal point flux */
214 
215  if((ithermal[0]>=2)&&(intpointvart==1)){
216 
217  NNEW(fn1,double,num_cpus*mt**nk);
218  NNEW(qa1,double,num_cpus*4);
219  NNEW(nal,ITG,num_cpus);
220 
221  co1=co;kon1=kon;ipkon1=ipkon;lakon1=lakon;v1=v;
222  elcon1=elcon;nelcon1=nelcon;rhcon1=rhcon;nrhcon1=nrhcon;
223  ielmat1=ielmat;ielorien1=ielorien;norien1=norien;orab1=orab;
224  ntmat1_=ntmat_;t01=t0;iperturb1=iperturb;iout1=iout;vold1=vold;
225  ipompc1=ipompc;nodempc1=nodempc;coefmpc1=coefmpc;nmpc1=nmpc;
226  dtime1=dtime;time1=time;ttime1=ttime;plkcon1=plkcon;
227  nplkcon1=nplkcon;xstateini1=xstateini;xstiff1=xstiff;
228  xstate1=xstate;npmat1_=npmat_;matname1=matname;mi1=mi;
229  ncmat1_=ncmat_;nstate1_=nstate_;cocon1=cocon;ncocon1=ncocon;
230  qfx1=qfx;ikmpc1=ikmpc;ilmpc1=ilmpc;istep1=istep;iinc1=iinc;
231  springarea1=springarea;calcul_fn1=calcul_fn;calcul_qa1=calcul_qa;
232  mt1=mt;nk1=nk;shcon1=shcon;nshcon1=nshcon;ithermal1=ithermal;
233  nelemload1=nelemload;nload1=nload;nmethod1=nmethod;reltime1=reltime;
234  sideload1=sideload;xload1=xload;xloadold1=xloadold;
235  pslavsurf1=pslavsurf;pmastsurf1=pmastsurf;mortar1=mortar;
236  clearini1=clearini;plicon1=plicon;nplicon1=nplicon;ielprop1=ielprop;
237  prop1=prop;iponoel1=iponoel;inoel1=inoel;network1=network;
238  ipobody1=ipobody;ibody1=ibody;xbody1=xbody;
239 
240  /* calculating the heat flux */
241 
242  printf(" Using up to %" ITGFORMAT " cpu(s) for the heat flux calculation.\n\n", num_cpus);
243 
244  /* create threads and wait */
245 
246  NNEW(ithread,ITG,num_cpus);
247  for(i=0; i<num_cpus; i++) {
248  ithread[i]=i;
249  pthread_create(&tid[i], NULL, (void *)resultsthermemmt, (void *)&ithread[i]);
250  }
251  for(i=0; i<num_cpus; i++) pthread_join(tid[i], NULL);
252 
253  for(i=0;i<*nk;i++){
254  fn[mt*i]=fn1[mt*i];
255  }
256  for(i=0;i<*nk;i++){
257  for(j=1;j<num_cpus;j++){
258  fn[mt*i]+=fn1[mt*i+j*mt**nk];
259  }
260  }
261  SFREE(fn1);SFREE(ithread);
262 
263  /* determine the internal concentrated heat flux */
264 
265  qa[1]=qa1[1];
266  for(j=1;j<num_cpus;j++){
267  qa[1]+=qa1[1+j*4];
268  }
269 
270  SFREE(qa1);
271 
272  for(j=1;j<num_cpus;j++){
273  nal[0]+=nal[j];
274  }
275 
276  if(calcul_qa==1){
277  if(nal[0]>0){
278  qa[1]/=nal[0];
279  }
280  }
281  SFREE(nal);
282  }
283 
284  /* calculating the thermal internal forces */
285 
286  FORTRAN(resultsforc_em,(nk,f,fn,nactdof,ipompc,nodempc,
287  coefmpc,labmpc,nmpc,mi,fmpc,&calcul_fn,&calcul_f,inomat));
288 
289  /* storing results in the .dat file
290  extrapolation of integration point values to the nodes
291  interpolation of 3d results for 1d/2d elements */
292 
293  FORTRAN(resultsprint,(co,nk,kon,ipkon,lakon,ne,v,stn,inum,
294  sti,ielorien,norien,orab,t1,ithermal,filab,een,iperturb,fn,
295  nactdof,iout,vold,nodeboun,ndirboun,nboun,nmethod,ttime,xstate,
296  epn,mi,
297  nstate_,ener,enern,xstaten,eei,set,nset,istartset,iendset,
298  ialset,nprint,prlab,prset,qfx,qfn,trab,inotr,ntrans,
299  nelemload,nload,&ikin,ielmat,thicke,eme,emn,rhcon,nrhcon,shcon,
300  nshcon,cocon,ncocon,ntmat_,sideload,icfd,inomat,pslavsurf,islavact,
301  cdn,&mortar,islavnode,nslavnode,ntie,islavsurf,time,ielprop,prop,
302  veold,ne0,nmpc,ipompc,nodempc,labmpc,energyini,energy,orname,
303  xload));
304 
305  return;
306 
307 }
static ITG * iinc1
Definition: resultsinduction.c:27
#define ITGFORMAT
Definition: CalculiX.h:52
static ITG * mi1
Definition: resultsinduction.c:27
static double * v1
Definition: resultsinduction.c:35
subroutine resultsini_em(nk, v, ithermal, filab, iperturb, f, fn, nactdof, iout, qa, b, nodeboun, ndirboun, xboun, nboun, ipompc, nodempc, coefmpc, labmpc, nmpc, nmethod, cam, neq, veold, dtime, mi, vini, nprint, prlab, intpointvarm, calcul_fn, calcul_f, calcul_qa, calcul_cauchy, nener, ikin, intpointvart, xforc, nforc)
Definition: resultsini_em.f:25
static ITG * nelemload1
Definition: resultsinduction.c:31
static ITG * nstate1_
Definition: resultsinduction.c:27
static double * prop1
Definition: resultsinduction.c:36
static double * qa1
Definition: resultsinduction.c:36
static ITG * ipobody1
Definition: resultsinduction.c:31
int pthread_create(pthread_t *thread_id, const pthread_attr_t *attributes, void *(*thread_function)(void *), void *arguments)
static double * sti1
Definition: resultsinduction.c:36
static ITG num_cpus
Definition: resultsinduction.c:31
static ITG calcul_fn1
Definition: resultsinduction.c:27
static ITG * istartset1
Definition: resultsinduction.c:31
static ITG * inoel1
Definition: resultsinduction.c:27
static double * elcon1
Definition: resultsinduction.c:35
static ITG * ielorien1
Definition: resultsinduction.c:27
static double * vini1
Definition: resultsinduction.c:36
static ITG * ibody1
Definition: resultsinduction.c:31
static double * xstate1
Definition: resultsinduction.c:36
static double * clearini1
Definition: resultsinduction.c:36
static char * sideload1
Definition: resultsinduction.c:25
static char * matname1
Definition: resultsinduction.c:25
static double * xloadold1
Definition: resultsinduction.c:36
static ITG * nelcon1
Definition: resultsinduction.c:27
static ITG * nload1
Definition: resultsinduction.c:31
static ITG * ielprop1
Definition: resultsinduction.c:27
void FORTRAN(actideacti,(char *set, ITG *nset, ITG *istartset, ITG *iendset, ITG *ialset, char *objectset, ITG *ipkon, ITG *ibject, ITG *ne))
static double * qfx1
Definition: resultsinduction.c:36
static ITG * kon1
Definition: resultsinduction.c:27
ITG strcmp1(const char *s1, const char *s2)
Definition: strcmp1.c:24
static ITG * iperturb1
Definition: resultsinduction.c:27
static double * alcon1
Definition: resultsinduction.c:35
static ITG mortar1
Definition: resultsinduction.c:31
subroutine frditeration(co, nk, kon, ipkon, lakon, ne, v, time, ielmat, matname, mi, istep, iinc, ithermal)
Definition: frditeration.f:21
static double * xbody1
Definition: resultsinduction.c:36
static ITG * nal
Definition: resultsinduction.c:31
ITG getSystemCPUs()
Definition: getSystemCPUs.c:40
static ITG * nodempc1
Definition: resultsinduction.c:31
static ITG * ipkon1
Definition: resultsinduction.c:27
static ITG * network1
Definition: resultsinduction.c:31
subroutine resultsprint(co, nk, kon, ipkon, lakon, ne, v, stn, inum, stx, ielorien, norien, orab, t1, ithermal, filab, een, iperturb, fn, nactdof, iout, vold, nodeboun, ndirboun, nboun, nmethod, ttime, xstate, epn, mi, nstate_, ener, enern, xstaten, eei, set, nset, istartset, iendset, ialset, nprint, prlab, prset, qfx, qfn, trab, inotr, ntrans, nelemload, nload, ikin, ielmat, thicke, eme, emn, rhcon, nrhcon, shcon, nshcon, cocon, ncocon, ntmat_, sideload, icfd, inomat, pslavsurf, islavact, cdn, mortar, islavnode, nslavnode, ntie, islavsurf, time, ielprop, prop, veold, ne0, nmpc, ipompc, nodempc, labmpc, energyini, energy, orname, xload)
Definition: resultsprint.f:29
static ITG * ipompc1
Definition: resultsinduction.c:31
static ITG calcul_qa1
Definition: resultsinduction.c:27
static double * springarea1
Definition: resultsinduction.c:36
static double * ttime1
Definition: resultsinduction.c:36
static ITG * nalcon1
Definition: resultsinduction.c:27
static ITG * iout1
Definition: resultsinduction.c:27
static ITG * iponoel1
Definition: resultsinduction.c:27
static ITG * nk1
Definition: resultsinduction.c:31
static ITG * nshcon1
Definition: resultsinduction.c:31
static ITG * nplkcon1
Definition: resultsinduction.c:27
static ITG * ithermal1
Definition: resultsinduction.c:27
static double * pmastsurf1
Definition: resultsinduction.c:36
#define SFREE(a)
Definition: CalculiX.h:41
static ITG * ncocon1
Definition: resultsinduction.c:31
static double * xstiff1
Definition: resultsinduction.c:36
static ITG * npmat1_
Definition: resultsinduction.c:27
static double * fn1
Definition: resultsinduction.c:36
static double * pslavsurf1
Definition: resultsinduction.c:36
static double * reltime1
Definition: resultsinduction.c:36
static ITG * nrhcon1
Definition: resultsinduction.c:27
static ITG * ielmat1
Definition: resultsinduction.c:27
subroutine resultsforc_em(nk, f, fn, nactdof, ipompc, nodempc, coefmpc, labmpc, nmpc, mi, fmpc, calcul_fn, calcul_f, inomat)
Definition: resultsforc_em.f:21
static double * t01
Definition: resultsinduction.c:35
static double * plicon1
Definition: resultsinduction.c:36
static ITG * iendset1
Definition: resultsinduction.c:31
static ITG * ikmpc1
Definition: resultsinduction.c:31
static double * xload1
Definition: resultsinduction.c:36
void * resultsthermemmt(ITG *i)
Definition: resultsinduction.c:330
static double * vold1
Definition: resultsinduction.c:36
static ITG * ilmpc1
Definition: resultsinduction.c:31
static ITG * istep1
Definition: resultsinduction.c:27
static char * lakon1
Definition: resultsinduction.c:25
int pthread_join(pthread_t thread, void **status_ptr)
static ITG * nmpc1
Definition: resultsinduction.c:31
static ITG * nmethod1
Definition: resultsinduction.c:27
static double * shcon1
Definition: resultsinduction.c:36
static double * orab1
Definition: resultsinduction.c:35
static ITG mt1
Definition: resultsinduction.c:31
#define ITG
Definition: CalculiX.h:51
static double * cocon1
Definition: resultsinduction.c:36
static double * time1
Definition: resultsinduction.c:36
static double * rhcon1
Definition: resultsinduction.c:35
static ITG * ialset1
Definition: resultsinduction.c:31
static ITG * ncmat1_
Definition: resultsinduction.c:27
static double * coefmpc1
Definition: resultsinduction.c:36
static double * co1
Definition: resultsinduction.c:35
void * resultsemmt(ITG *i)
Definition: resultsinduction.c:311
static ITG * nplicon1
Definition: resultsinduction.c:27
static double * xstateini1
Definition: resultsinduction.c:36
#define NNEW(a, b, c)
Definition: CalculiX.h:39
static ITG * ne1
Definition: resultsinduction.c:27
static double * dtime1
Definition: resultsinduction.c:36
static ITG * ntmat1_
Definition: resultsinduction.c:27
static ITG * iactive1
Definition: resultsinduction.c:31
static ITG * norien1
Definition: resultsinduction.c:27
static double * plkcon1
Definition: resultsinduction.c:36
static double * h01
Definition: resultsinduction.c:36

◆ resultsthermemmt()

void* resultsthermemmt ( ITG i)
330  {
331 
332  ITG indexfn,indexqa,indexnal,nea,neb,nedelta;
333 
334  indexfn=*i*mt1**nk1;
335  indexqa=*i*4;
336  indexnal=*i;
337 
338  nedelta=(ITG)floor(*ne1/(double)num_cpus);
339  nea=*i*nedelta+1;
340  neb=(*i+1)*nedelta;
341  if((*i==num_cpus-1)&&(neb<*ne1)) neb=*ne1;
342 
348  npmat1_,
351  &calcul_fn1,&calcul_qa1,&nal[indexnal],&nea,&neb,ithermal1,
355 
356  return NULL;
357 }
static ITG * iinc1
Definition: resultsinduction.c:27
static ITG * mi1
Definition: resultsinduction.c:27
static double * v1
Definition: resultsinduction.c:35
static ITG * nelemload1
Definition: resultsinduction.c:31
static ITG * nstate1_
Definition: resultsinduction.c:27
static double * prop1
Definition: resultsinduction.c:36
static double * qa1
Definition: resultsinduction.c:36
static ITG * ipobody1
Definition: resultsinduction.c:31
static ITG num_cpus
Definition: resultsinduction.c:31
static ITG calcul_fn1
Definition: resultsinduction.c:27
static ITG * inoel1
Definition: resultsinduction.c:27
static double * elcon1
Definition: resultsinduction.c:35
static ITG * ielorien1
Definition: resultsinduction.c:27
static ITG * ibody1
Definition: resultsinduction.c:31
static double * xstate1
Definition: resultsinduction.c:36
static double * clearini1
Definition: resultsinduction.c:36
static char * sideload1
Definition: resultsinduction.c:25
subroutine resultstherm(co, kon, ipkon, lakon, v, elcon, nelcon, rhcon, nrhcon, ielmat, ielorien, norien, orab, ntmat_, t0, iperturb, fn, shcon, nshcon, iout, qa, vold, ipompc, nodempc, coefmpc, nmpc, dtime, time, ttime, plkcon, nplkcon, xstateini, xstiff, xstate, npmat_, matname, mi, ncmat_, nstate_, cocon, ncocon, qfx, ikmpc, ilmpc, istep, iinc, springarea, calcul_fn, calcul_qa, nal, nea, neb, ithermal, nelemload, nload, nmethod, reltime, sideload, xload, xloadold, pslavsurf, pmastsurf, mortar, clearini, plicon, nplicon, ielprop, prop, iponoel, inoel, network, ipobody, xbody, ibody)
Definition: resultstherm.f:30
static char * matname1
Definition: resultsinduction.c:25
static double * xloadold1
Definition: resultsinduction.c:36
static ITG * nelcon1
Definition: resultsinduction.c:27
static ITG * nload1
Definition: resultsinduction.c:31
static ITG * ielprop1
Definition: resultsinduction.c:27
void FORTRAN(actideacti,(char *set, ITG *nset, ITG *istartset, ITG *iendset, ITG *ialset, char *objectset, ITG *ipkon, ITG *ibject, ITG *ne))
static double * qfx1
Definition: resultsinduction.c:36
static ITG * kon1
Definition: resultsinduction.c:27
static ITG * iperturb1
Definition: resultsinduction.c:27
static ITG mortar1
Definition: resultsinduction.c:31
static double * xbody1
Definition: resultsinduction.c:36
static ITG * nal
Definition: resultsinduction.c:31
static ITG * nodempc1
Definition: resultsinduction.c:31
static ITG * ipkon1
Definition: resultsinduction.c:27
static ITG * network1
Definition: resultsinduction.c:31
static ITG * ipompc1
Definition: resultsinduction.c:31
static ITG calcul_qa1
Definition: resultsinduction.c:27
static double * springarea1
Definition: resultsinduction.c:36
static double * ttime1
Definition: resultsinduction.c:36
static ITG * iout1
Definition: resultsinduction.c:27
static ITG * iponoel1
Definition: resultsinduction.c:27
static ITG * nk1
Definition: resultsinduction.c:31
static ITG * nshcon1
Definition: resultsinduction.c:31
static ITG * nplkcon1
Definition: resultsinduction.c:27
static ITG * ithermal1
Definition: resultsinduction.c:27
static double * pmastsurf1
Definition: resultsinduction.c:36
static ITG * ncocon1
Definition: resultsinduction.c:31
static double * xstiff1
Definition: resultsinduction.c:36
static ITG * npmat1_
Definition: resultsinduction.c:27
static double * fn1
Definition: resultsinduction.c:36
static double * pslavsurf1
Definition: resultsinduction.c:36
static double * reltime1
Definition: resultsinduction.c:36
static ITG * nrhcon1
Definition: resultsinduction.c:27
static ITG * ielmat1
Definition: resultsinduction.c:27
static double * t01
Definition: resultsinduction.c:35
static double * plicon1
Definition: resultsinduction.c:36
static ITG * ikmpc1
Definition: resultsinduction.c:31
static double * xload1
Definition: resultsinduction.c:36
static double * vold1
Definition: resultsinduction.c:36
static ITG * ilmpc1
Definition: resultsinduction.c:31
static ITG * istep1
Definition: resultsinduction.c:27
static char * lakon1
Definition: resultsinduction.c:25
static ITG * nmpc1
Definition: resultsinduction.c:31
static ITG * nmethod1
Definition: resultsinduction.c:27
static double * shcon1
Definition: resultsinduction.c:36
static double * orab1
Definition: resultsinduction.c:35
static ITG mt1
Definition: resultsinduction.c:31
#define ITG
Definition: CalculiX.h:51
static double * cocon1
Definition: resultsinduction.c:36
static double * time1
Definition: resultsinduction.c:36
static double * rhcon1
Definition: resultsinduction.c:35
static ITG * ncmat1_
Definition: resultsinduction.c:27
static double * coefmpc1
Definition: resultsinduction.c:36
static double * co1
Definition: resultsinduction.c:35
static ITG * nplicon1
Definition: resultsinduction.c:27
static double * xstateini1
Definition: resultsinduction.c:36
static ITG * ne1
Definition: resultsinduction.c:27
static double * dtime1
Definition: resultsinduction.c:36
static ITG * ntmat1_
Definition: resultsinduction.c:27
static ITG * norien1
Definition: resultsinduction.c:27
static double * plkcon1
Definition: resultsinduction.c:36

Variable Documentation

◆ alcon1

double * alcon1
static

◆ calcul_fn1

ITG calcul_fn1
static

◆ calcul_qa1

ITG calcul_qa1
static

◆ clearini1

double * clearini1
static

◆ co1

double* co1
static

◆ cocon1

double * cocon1
static

◆ coefmpc1

double * coefmpc1
static

◆ dtime1

double * dtime1
static

◆ elcon1

double * elcon1
static

◆ fn1

double * fn1 =NULL
static

◆ h01

double * h01
static

◆ iactive1

ITG * iactive1
static

◆ ialset1

ITG * ialset1
static

◆ ibody1

ITG * ibody1
static

◆ ielmat1

ITG * ielmat1
static

◆ ielorien1

ITG * ielorien1
static

◆ ielprop1

ITG * ielprop1
static

◆ iendset1

ITG * iendset1
static

◆ iinc1

ITG * iinc1
static

◆ ikmpc1

ITG * ikmpc1
static

◆ ilmpc1

ITG * ilmpc1
static

◆ inoel1

ITG * inoel1
static

◆ iout1

ITG * iout1
static

◆ iperturb1

ITG * iperturb1
static

◆ ipkon1

ITG * ipkon1
static

◆ ipobody1

ITG * ipobody1
static

◆ ipompc1

ITG * ipompc1
static

◆ iponoel1

ITG * iponoel1
static

◆ istartset1

ITG * istartset1
static

◆ istep1

ITG * istep1
static

◆ ithermal1

ITG * ithermal1
static

◆ kon1

ITG* kon1
static

◆ lakon1

char* lakon1
static

◆ matname1

char * matname1
static

◆ mi1

ITG * mi1
static

◆ mortar1

ITG mortar1
static

◆ mt1

ITG mt1
static

◆ nal

ITG * nal =NULL
static

◆ nalcon1

ITG * nalcon1
static

◆ ncmat1_

ITG * ncmat1_
static

◆ ncocon1

ITG * ncocon1
static

◆ ne1

ITG * ne1
static

◆ nelcon1

ITG * nelcon1
static

◆ nelemload1

ITG * nelemload1
static

◆ network1

ITG * network1
static

◆ nk1

ITG * nk1
static

◆ nload1

ITG * nload1
static

◆ nmethod1

ITG * nmethod1
static

◆ nmpc1

ITG * nmpc1
static

◆ nodempc1

ITG * nodempc1
static

◆ norien1

ITG * norien1
static

◆ nplicon1

ITG * nplicon1
static

◆ nplkcon1

ITG * nplkcon1
static

◆ npmat1_

ITG * npmat1_
static

◆ nrhcon1

ITG * nrhcon1
static

◆ nshcon1

ITG * nshcon1
static

◆ nstate1_

ITG * nstate1_
static

◆ ntmat1_

ITG * ntmat1_
static

◆ num_cpus

ITG num_cpus
static

◆ orab1

double * orab1
static

◆ plicon1

double * plicon1
static

◆ plkcon1

double * plkcon1
static

◆ pmastsurf1

double * pmastsurf1
static

◆ prop1

double * prop1
static

◆ pslavsurf1

double * pslavsurf1
static

◆ qa1

double * qa1 =NULL
static

◆ qfx1

double * qfx1
static

◆ reltime1

double * reltime1
static

◆ rhcon1

double * rhcon1
static

◆ shcon1

double * shcon1
static

◆ sideload1

char * sideload1
static

◆ springarea1

double * springarea1
static

◆ sti1

double * sti1
static

◆ t01

double * t01
static

◆ time1

double * time1
static

◆ ttime1

double * ttime1
static

◆ v1

double * v1
static

◆ vini1

double * vini1
static

◆ vold1

double * vold1
static

◆ xbody1

double * xbody1
static

◆ xload1

double * xload1
static

◆ xloadold1

double * xloadold1
static

◆ xstate1

double * xstate1
static

◆ xstateini1

double * xstateini1
static

◆ xstiff1

double * xstiff1
static
Hosted by OpenAircraft.com, (Michigan UAV, LLC)