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

Go to the source code of this file.

Functions

void mafillvmain (ITG *nef, ITG *ipnei, ITG *neifa, ITG *neiel, double *vfa, double *xxn, double *area, double *auv, double *adv, ITG *jq, ITG *irow, ITG *nzs, double *bv, double *vel, double *cosa, double *umfa, double *xlet, double *xle, double *gradvfa, double *xxi, double *body, double *volume, ITG *ielfa, char *lakonf, ITG *ifabou, ITG *nbody, double *dtimef, double *velo, double *veloo, double *sel, double *xrlfa, double *gamma, double *xxj, ITG *nactdohinv, double *a1, double *a2, double *a3, double *flux, ITG *icyclic, double *c, ITG *ifatie, ITG *iau6, double *xxni, double *xxnj, ITG *iturbulent, double *gradvel)
 
void * mafillvmt (ITG *i)
 

Variables

static char * lakonf1
 
static ITG num_cpus
 
static ITGnef1
 
static ITGipnei1
 
static ITGneifa1
 
static ITGneiel1
 
static ITGjq1
 
static ITGirow1
 
static ITGnzs1
 
static ITGielfa1
 
static ITGifabou1
 
static ITGnbody1
 
static ITGnactdohinv1
 
static ITGicyclic1
 
static ITGifatie1
 
static ITGiau61
 
static ITGiturbulent1
 
static double * auv1
 
static double * adv1
 
static double * bv1
 
static double * vfa1
 
static double * xxn1
 
static double * area1
 
static double * vel1
 
static double * cosa1
 
static double * umfa1
 
static double * xlet1
 
static double * xle1
 
static double * gradvfa1
 
static double * xxi1
 
static double * body1
 
static double * volume1
 
static double * dtimef1
 
static double * velo1
 
static double * veloo1
 
static double * sel1
 
static double * xrlfa1
 
static double * gamma1
 
static double * xxj1
 
static double * a11
 
static double * a21
 
static double * a31
 
static double * flux1
 
static double * c1
 
static double * xxni1
 
static double * xxnj1
 
static double * gradvel1
 

Function Documentation

◆ mafillvmain()

void mafillvmain ( ITG nef,
ITG ipnei,
ITG neifa,
ITG neiel,
double *  vfa,
double *  xxn,
double *  area,
double *  auv,
double *  adv,
ITG jq,
ITG irow,
ITG nzs,
double *  bv,
double *  vel,
double *  cosa,
double *  umfa,
double *  xlet,
double *  xle,
double *  gradvfa,
double *  xxi,
double *  body,
double *  volume,
ITG ielfa,
char *  lakonf,
ITG ifabou,
ITG nbody,
double *  dtimef,
double *  velo,
double *  veloo,
double *  sel,
double *  xrlfa,
double *  gamma,
double *  xxj,
ITG nactdohinv,
double *  a1,
double *  a2,
double *  a3,
double *  flux,
ITG icyclic,
double *  c,
ITG ifatie,
ITG iau6,
double *  xxni,
double *  xxnj,
ITG iturbulent,
double *  gradvel 
)
45  {
46 
47  ITG i,j;
48 
49  /* variables for multithreading procedure */
50 
51  ITG sys_cpus,*ithread=NULL;
52  char *env,*envloc,*envsys;
53 
54  num_cpus = 0;
55  sys_cpus=0;
56 
57  /* explicit user declaration prevails */
58 
59  envsys=getenv("NUMBER_OF_CPUS");
60  if(envsys){
61  sys_cpus=atoi(envsys);
62  if(sys_cpus<0) sys_cpus=0;
63  }
64 
65  /* automatic detection of available number of processors */
66 
67  if(sys_cpus==0){
68  sys_cpus = getSystemCPUs();
69  if(sys_cpus<1) sys_cpus=1;
70  }
71 
72  /* local declaration prevails, if strictly positive */
73 
74  envloc = getenv("CCX_NPROC_CFD");
75  if(envloc){
76  num_cpus=atoi(envloc);
77  if(num_cpus<0){
78  num_cpus=0;
79  }else if(num_cpus>sys_cpus){
80  num_cpus=sys_cpus;
81  }
82 
83  }
84 
85  /* else global declaration, if any, applies */
86 
87  env = getenv("OMP_NUM_THREADS");
88  if(num_cpus==0){
89  if (env)
90  num_cpus = atoi(env);
91  if (num_cpus < 1) {
92  num_cpus=1;
93  }else if(num_cpus>sys_cpus){
94  num_cpus=sys_cpus;
95  }
96  }
97 
98 // next line is to be inserted in a similar way for all other paralell parts
99 
100  if(*nef<num_cpus) num_cpus=*nef;
101 
102  pthread_t tid[num_cpus];
103 
104  /* calculating the stiffness and/or mass matrix
105  (symmetric part) */
106 
107  nef1=nef;ipnei1=ipnei;neifa1=neifa;neiel1=neiel;vfa1=vfa;xxn1=xxn;
108  area1=area;jq1=jq;irow1=irow;nzs1=nzs;vel1=vel;cosa1=cosa;umfa1=umfa;
109  xlet1=xlet;xle1=xle;gradvfa1=gradvfa;xxi1=xxi;body1=body;volume1=volume;
110  ielfa1=ielfa;lakonf1=lakonf;ifabou1=ifabou;nbody1=nbody;
111  dtimef1=dtimef;velo1=velo;veloo1=veloo;sel1=sel;xrlfa1=xrlfa;
112  gamma1=gamma;xxj1=xxj;nactdohinv1=nactdohinv;a11=a1;a21=a2;a31=a3;
113  flux1=flux;icyclic1=icyclic;c1=c;ifatie1=ifatie;iau61=iau6;
114  adv1=adv;auv1=auv;bv1=bv;xxni1=xxni;xxnj1=xxnj;iturbulent1=iturbulent;
115  gradvel1=gradvel;
116 
117  /* create threads and wait */
118 
119  NNEW(ithread,ITG,num_cpus);
120  for(i=0; i<num_cpus; i++) {
121  ithread[i]=i;
122  pthread_create(&tid[i], NULL, (void *)mafillvmt, (void *)&ithread[i]);
123  }
124  for(i=0; i<num_cpus; i++) pthread_join(tid[i], NULL);
125 
126  SFREE(ithread);
127 
128  return;
129 
130 }
static double * volume1
Definition: mafillvmain.c:30
static double * gradvel1
Definition: mafillvmain.c:30
static ITG * ipnei1
Definition: mafillvmain.c:27
static double * xle1
Definition: mafillvmain.c:30
static ITG * irow1
Definition: mafillvmain.c:27
void * mafillvmt(ITG *i)
Definition: mafillvmain.c:134
int pthread_create(pthread_t *thread_id, const pthread_attr_t *attributes, void *(*thread_function)(void *), void *arguments)
static double * dtimef1
Definition: mafillvmain.c:30
static ITG * icyclic1
Definition: mafillvmain.c:27
static double * auv1
Definition: mafillvmain.c:30
static double * cosa1
Definition: mafillvmain.c:30
subroutine flux(node1, node2, nodem, nelem, lakon, kon, ipkon, nactdog, identity, ielprop, prop, kflag, v, xflow, f, nodef, idirf, df, cp, R, rho, physcon, g, co, dvi, numf, vold, set, shcon, nshcon, rhcon, nrhcon, ntmat_, mi, ider, ttime, time, iaxial)
Definition: flux.f:24
static ITG num_cpus
Definition: mafillvmain.c:27
static double * area1
Definition: mafillvmain.c:30
static ITG * nzs1
Definition: mafillvmain.c:27
static ITG * ifatie1
Definition: mafillvmain.c:27
static double * xxni1
Definition: mafillvmain.c:30
ITG getSystemCPUs()
Definition: getSystemCPUs.c:40
static double * a21
Definition: mafillvmain.c:30
static double * vfa1
Definition: mafillvmain.c:30
static double * a11
Definition: mafillvmain.c:30
static double * flux1
Definition: mafillvmain.c:30
static ITG * nbody1
Definition: mafillvmain.c:27
static ITG * iau61
Definition: mafillvmain.c:27
static double * xxj1
Definition: mafillvmain.c:30
static double * velo1
Definition: mafillvmain.c:30
static double * body1
Definition: mafillvmain.c:30
static double * gradvfa1
Definition: mafillvmain.c:30
static ITG * ielfa1
Definition: mafillvmain.c:27
static double * umfa1
Definition: mafillvmain.c:30
static double * bv1
Definition: mafillvmain.c:30
#define SFREE(a)
Definition: CalculiX.h:41
static ITG * ifabou1
Definition: mafillvmain.c:27
static double * a31
Definition: mafillvmain.c:30
static ITG * jq1
Definition: mafillvmain.c:27
static ITG * nactdohinv1
Definition: mafillvmain.c:27
static double * sel1
Definition: mafillvmain.c:30
static double * xxnj1
Definition: mafillvmain.c:30
static double * xrlfa1
Definition: mafillvmain.c:30
static double * vel1
Definition: mafillvmain.c:30
static double * xxi1
Definition: mafillvmain.c:30
static double * c1
Definition: mafillvmain.c:30
static double * gamma1
Definition: mafillvmain.c:30
int pthread_join(pthread_t thread, void **status_ptr)
static double * xlet1
Definition: mafillvmain.c:30
static double * veloo1
Definition: mafillvmain.c:30
static ITG * neiel1
Definition: mafillvmain.c:27
static char * lakonf1
Definition: mafillvmain.c:25
#define ITG
Definition: CalculiX.h:51
static double * adv1
Definition: mafillvmain.c:30
#define NNEW(a, b, c)
Definition: CalculiX.h:39
static double * xxn1
Definition: mafillvmain.c:30
static ITG * neifa1
Definition: mafillvmain.c:27
static ITG * nef1
Definition: mafillvmain.c:27
static ITG * iturbulent1
Definition: mafillvmain.c:27

◆ mafillvmt()

void* mafillvmt ( ITG i)
134  {
135 
136  ITG nefa,nefb,nefdelta;
137 
138 // ceil -> floor
139 
140  nefdelta=(ITG)floor(*nef1/(double)num_cpus);
141  nefa=*i*nefdelta+1;
142  nefb=(*i+1)*nefdelta;
143 // next line! -> all parallel sections
144  if((*i==num_cpus-1)&&(nefb<*nef1)) nefb=*nef1;
145 
151  a11,a21,a31,flux1,&nefa,&nefb,icyclic1,c1,ifatie1,iau61,
153 
154  return NULL;
155 }
static double * volume1
Definition: mafillvmain.c:30
static double * gradvel1
Definition: mafillvmain.c:30
static ITG * ipnei1
Definition: mafillvmain.c:27
static double * xle1
Definition: mafillvmain.c:30
static ITG * irow1
Definition: mafillvmain.c:27
static double * dtimef1
Definition: mafillvmain.c:30
static ITG * icyclic1
Definition: mafillvmain.c:27
subroutine mafillv(nef, ipnei, neifa, neiel, vfa, xxn, area, auv, adv, jq, irow, nzs, bv, vel, cosa, umfa, xlet, xle, gradvfa, xxi, body, volume, ielfa, lakonf, ifabou, nbody, dtimef, velo, veloo, sel, xrlfa, gamma, xxj, nactdohinv, a1, a2, a3, flux, nefa, nefb, icyclic, c, ifatie, iau6, xxni, xxnj, iturbulent, gradvel)
Definition: mafillv.f:25
static double * auv1
Definition: mafillvmain.c:30
static double * cosa1
Definition: mafillvmain.c:30
void FORTRAN(actideacti,(char *set, ITG *nset, ITG *istartset, ITG *iendset, ITG *ialset, char *objectset, ITG *ipkon, ITG *ibject, ITG *ne))
static ITG num_cpus
Definition: mafillvmain.c:27
static double * area1
Definition: mafillvmain.c:30
static ITG * nzs1
Definition: mafillvmain.c:27
static ITG * ifatie1
Definition: mafillvmain.c:27
static double * xxni1
Definition: mafillvmain.c:30
static double * a21
Definition: mafillvmain.c:30
static double * vfa1
Definition: mafillvmain.c:30
static double * a11
Definition: mafillvmain.c:30
static double * flux1
Definition: mafillvmain.c:30
static ITG * nbody1
Definition: mafillvmain.c:27
static ITG * iau61
Definition: mafillvmain.c:27
static double * xxj1
Definition: mafillvmain.c:30
static double * velo1
Definition: mafillvmain.c:30
static double * body1
Definition: mafillvmain.c:30
static double * gradvfa1
Definition: mafillvmain.c:30
static ITG * ielfa1
Definition: mafillvmain.c:27
static double * umfa1
Definition: mafillvmain.c:30
static double * bv1
Definition: mafillvmain.c:30
static ITG * ifabou1
Definition: mafillvmain.c:27
static double * a31
Definition: mafillvmain.c:30
static ITG * jq1
Definition: mafillvmain.c:27
static ITG * nactdohinv1
Definition: mafillvmain.c:27
static double * sel1
Definition: mafillvmain.c:30
static double * xxnj1
Definition: mafillvmain.c:30
static double * xrlfa1
Definition: mafillvmain.c:30
static double * vel1
Definition: mafillvmain.c:30
static double * xxi1
Definition: mafillvmain.c:30
static double * c1
Definition: mafillvmain.c:30
static double * gamma1
Definition: mafillvmain.c:30
static double * xlet1
Definition: mafillvmain.c:30
static double * veloo1
Definition: mafillvmain.c:30
static ITG * neiel1
Definition: mafillvmain.c:27
static char * lakonf1
Definition: mafillvmain.c:25
#define ITG
Definition: CalculiX.h:51
static double * adv1
Definition: mafillvmain.c:30
static double * xxn1
Definition: mafillvmain.c:30
static ITG * neifa1
Definition: mafillvmain.c:27
static ITG * nef1
Definition: mafillvmain.c:27
static ITG * iturbulent1
Definition: mafillvmain.c:27

Variable Documentation

◆ a11

double * a11
static

◆ a21

double * a21
static

◆ a31

double * a31
static

◆ adv1

double * adv1
static

◆ area1

double * area1
static

◆ auv1

double* auv1
static

◆ body1

double * body1
static

◆ bv1

double * bv1
static

◆ c1

double * c1
static

◆ cosa1

double * cosa1
static

◆ dtimef1

double * dtimef1
static

◆ flux1

double * flux1
static

◆ gamma1

double * gamma1
static

◆ gradvel1

double * gradvel1
static

◆ gradvfa1

double * gradvfa1
static

◆ iau61

ITG * iau61
static

◆ icyclic1

ITG * icyclic1
static

◆ ielfa1

ITG * ielfa1
static

◆ ifabou1

ITG * ifabou1
static

◆ ifatie1

ITG * ifatie1
static

◆ ipnei1

ITG * ipnei1
static

◆ irow1

ITG * irow1
static

◆ iturbulent1

ITG * iturbulent1
static

◆ jq1

ITG * jq1
static

◆ lakonf1

char* lakonf1
static

◆ nactdohinv1

ITG * nactdohinv1
static

◆ nbody1

ITG * nbody1
static

◆ nef1

ITG * nef1
static

◆ neiel1

ITG * neiel1
static

◆ neifa1

ITG * neifa1
static

◆ num_cpus

ITG num_cpus
static

◆ nzs1

ITG * nzs1
static

◆ sel1

double * sel1
static

◆ umfa1

double * umfa1
static

◆ vel1

double * vel1
static

◆ velo1

double * velo1
static

◆ veloo1

double * veloo1
static

◆ vfa1

double * vfa1
static

◆ volume1

double * volume1
static

◆ xle1

double * xle1
static

◆ xlet1

double * xlet1
static

◆ xrlfa1

double * xrlfa1
static

◆ xxi1

double * xxi1
static

◆ xxj1

double * xxj1
static

◆ xxn1

double * xxn1
static

◆ xxni1

double * xxni1
static

◆ xxnj1

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