AnyType yams_Op::operator () (Stack stack) const { // initialisation MeshPoint *mp(MeshPointStack(stack)), mps = *mp; Mesh3 *pTh = GetAny<Mesh3 *>((*eTh)(stack)); ffassert(pTh); Mesh3 &Th3 = *pTh; int nv = Th3.nv; int nt = Th3.nt; int nbe = Th3.nbe; KN<int> defaultintopt(23); KN<double> defaultfopt(14); defaultintopt = 0; defaultfopt = 0.; yams_inival(defaultintopt, defaultfopt); KN<int> intopt(23); for (int ii = 0; ii < 23; ii++) { intopt[ii] = defaultintopt[ii]; } KN<double> fopt(14); for (int ii = 0; ii < 14; ii++) { fopt[ii] = defaultfopt[ii]; } assert(fopt.N() == 14); if (nargs[0]) { KN<int> intopttmp = GetAny<KN_<long> >((*nargs[0])(stack)); if (intopttmp.N() != 13) { cerr << "the size of vector loptions is 13 " << endl; exit(1); } else { for (int ii = 0; ii < 13; ii++) { intopt[wrapper_intopt[ii]] = intopttmp[ii]; } } } if (nargs[1]) { KN<double> fopttmp = GetAny<KN_<double> >((*nargs[1])(stack)); if (fopttmp.N() != 11) { cerr << "the size of vector loptions is 11 not " << fopttmp.N() << endl; ExecError("FreeYams"); } else { for (int ii = 0; ii < 11; ii++) { fopt[wrapper_fopt[ii]] = fopttmp[ii]; } } } intopt[0] = arg(3, stack, intopt[0] != 1); intopt[8] = arg(4, stack, intopt[8]); fopt[7] = arg(5, stack, fopt[7]); fopt[8] = arg(6, stack, fopt[7]); fopt[6] = arg(7, stack, fopt[6]); intopt[22] = arg(8, stack, intopt[22]); // optim option if (nargs[9]) {intopt[17] = 1;} fopt[13] = arg(9, stack, fopt[13]); // ridge angle intopt[21] = arg(10, stack, intopt[21]);// absolue intopt[11] = arg(11, stack, (int)verbosity);// verbosity intopt[17] = arg(12, stack, intopt[17]);// no ridge intopt[18] = arg(13, stack, intopt[18]);// nb smooth if (verbosity > 1) { cout << " fopt = ["; for (int i = 0; i < 11; ++i) { cout << fopt[wrapper_fopt[i]] << (i < 10 ? "," : "];\n"); } cout << " intopt = ["; for (int i = 0; i < 13; ++i) { cout << intopt[wrapper_intopt[i]] << (i < 12 ? "," : "];\n"); } } /* * KN<int> intopt(arg(0,stack,defaultintopt)); * assert( intopt.N() == 23 ); * KN<double> fopt(arg(1,stack,defaultfopt)); * assert( fopt.N() == 14 ); */ KN<double> metric; int mtype = type; if (nargs[2]) { metric = GetAny<KN_<double> >((*nargs[2])(stack)); if (metric.N() == Th3.nv) { mtype = 1; intopt[1] = 0; } else if (metric.N() == 6 * Th3.nv) { intopt[1] = 1; mtype = 3; } else { cerr << "sizeof vector metric is incorrect, size will be Th.nv or 6*Th.nv" << endl; } } else if (nbsol > 0) { if (type == 1) { intopt[1] = 0; metric.resize(Th3.nv); metric = 0.; } else if (type == 3) { intopt[1] = 1; metric.resize(6 * Th3.nv); metric = 0.; } } else { if (intopt[1] == 0) {metric.resize(Th3.nv); metric = 0.;} else if (intopt[1] == 1) {metric.resize(6 * Th3.nv); metric = 0.;} } // mesh for yams yams_pSurfMesh yamsmesh; yamsmesh = (yams_pSurfMesh)calloc(1, sizeof(yams_SurfMesh)); if (!yamsmesh) { cerr << "allocation error for SurfMesh for yams" << endl; } yamsmesh->infile = NULL; yamsmesh->outfile = NULL; yamsmesh->type = M_SMOOTH | M_QUERY | M_DETECT | M_BINARY | M_OUTPUT; mesh3_to_yams_pSurfMesh(Th3, intopt[8], intopt[22], yamsmesh); // solution for freeyams2 if (nbsol) { MeshPoint *mp3(MeshPointStack(stack)); KN<bool> takemesh(nv); takemesh = false; for (int it = 0; it < nt; it++) { for (int iv = 0; iv < 4; iv++) { int i = Th3(it, iv); if (takemesh[i] == false) { mp3->setP(&Th3, it, iv); for (int ii = 0; ii < nbsolsize; ii++) { metric[i * nbsolsize + ii] = GetAny<double>((*sol[ii])(stack)); } takemesh[i] = true; } } } } if (verbosity > 10) { cout << "nbsol " << nargs[2] << endl; } if (nargs[2] || (nbsol > 0)) { float hmin, hmax; solyams_pSurfMesh(yamsmesh, mtype, metric, hmin, hmax); yamsmesh->nmfixe = yamsmesh->npfixe; if (fopt[7] < 0.0) { fopt[7] = max(fopt[7], hmin); } if (fopt[8] < 0.0) { fopt[8] = max(fopt[8], hmax); } } else { yamsmesh->nmfixe = 0; } int infondang = 0, infocc = 0; int res = yams_main(yamsmesh, intopt, fopt, infondang, infocc); if (verbosity > 10) { cout << " yamsmesh->dim " << yamsmesh->dim << endl; } if (res > 0) { cout << " problem with yams :: error " << res << endl; ExecError("Freeyams error"); } Mesh3 *Th3_T = yams_pSurfMesh_to_mesh3(yamsmesh, infondang, infocc, intopt[22]); // recuperer la solution ???? if (verbosity > 10) { cout << &yamsmesh->point << " " << &yamsmesh->tria << " " << &yamsmesh->geom << " " << &yamsmesh->tgte << endl; cout << &yamsmesh << endl; } free(yamsmesh->point); free(yamsmesh->tria); free(yamsmesh->geom); free(yamsmesh->tgte); if (yamsmesh->metric) {free(yamsmesh->metric);} if (yamsmesh->edge) {free(yamsmesh->edge);} if (yamsmesh->tetra) {free(yamsmesh->tetra);} free(yamsmesh); *mp = mps; Add2StackOfPtr2FreeRC(stack, Th3_T); return SetAny<pmesh3>(Th3_T); }
dSolvepastixmpi(const MatriceMorse<double> &AA, string datafile, KN<long> ¶m_int, KN<double> ¶m_double, KN<long> &pperm_r, KN<long> &pperm_c) : data_option(datafile) { //int m; //int ierr; struct timeval tv1, tv2; ia = NULL; ja = NULL; avals = NULL; loc2glob = NULL; rhs = NULL; pastix_data = NULL; // matrix assembled on host MPI_Comm_rank(MPI_COMM_WORLD, &myid); printf("- Rang MPI : %d\n", myid); MPI_Comm_size(MPI_COMM_WORLD, &mpi_size); // SYMETRIQUE mpi_flag = 0; thrd_flag = 0; Ncol = AA.m; Nrow = AA.n; // Avant : on ecrit la transposée /* ia = (pastix_int_t *) malloc( (Ncol+1)*sizeof(pastix_int_t)); for(int ii=0; ii < Ncol+1; ii++){ ia[ii] = AA.lg[ii]+1; } assert( ia[Ncol]-1 == AA.nbcoef ); ja = (pastix_int_t *) malloc((ia[Ncol]-1)*sizeof(pastix_int_t)); for(int ii=0; ii < ia[Ncol]-1; ii++) ja[ii] = AA.cl[ii]+1; if( sizeof(pastix_float_t) == sizeof(double) ){ avals = (pastix_float_t *) malloc( (ia[Ncol]-1)*sizeof(pastix_float_t)); for(int ii=0; ii < ia[Ncol]-1; ii++) avals[ii] = AA.a[ii]; } */ // AA.cl : indices des colonnes // AA.lg : pointeurs des lignes Morse_to_CSC( AA.n , AA.m, AA.nbcoef, AA.a, AA.cl, AA.lg, &avals, &ja, &ia); // ia : pointeurs des colonnes // ja : indices des lignes cout << "AA.n= "<< AA.n << " AA.m=" << AA.m << " AA.nbcoef=" << AA.nbcoef << endl; for(int ii=0; ii < Ncol+1; ii++){ ia[ii] = ia[ii]+1; } assert( ia[Ncol]-1 == AA.nbcoef ); for(int ii=0; ii < ia[Ncol]-1; ii++){ ja[ii] = ja[ii]+1; } perm = (pastix_int_t *) malloc(Ncol*sizeof(pastix_int_t)); invp = (pastix_int_t *) malloc(Ncol*sizeof(pastix_int_t)); rhs = (pastix_float_t *) malloc(Ncol*sizeof(pastix_float_t)); // reading permutation given by the user if(pperm_r) for(int ii=0; ii < Ncol; ii++) perm[ii] = pperm_r[ii]; if(pperm_c) for(int ii=0; ii < Ncol; ii++) invp[ii] = pperm_c[ii]; // CAS DE LA MATRICE NON DISTRIBUER pastix_int_t init_raff; fprintf(stdout,"-- INIT PARAMETERS --\n"); // reading iparm from array if(!data_option.empty()) read_datafile_pastixff(data_option,iparm,dparm); else if(param_int || param_double){ if( param_int ) { cout << "read param_int" << endl; assert(param_int.N() == 64); for(int ii=0; ii<64; ii++) iparm[ii] = param_int[ii]; iparm[IPARM_MODIFY_PARAMETER] = API_YES; } if( param_double ) { cout << "read param_double" << endl; assert(param_double.N() == 64); for(int ii=0; ii<64; ii++) dparm[ii] = param_double[ii]; } } else{ iparm[IPARM_MODIFY_PARAMETER] = API_NO; cout << "initialize parameter" << endl; } iparm[IPARM_START_TASK] = API_TASK_INIT; iparm[IPARM_END_TASK] = API_TASK_INIT; iparm[IPARM_SYM] = API_SYM_NO; // Matrix is considered nonsymetric pastix(&pastix_data, MPI_COMM_WORLD, Ncol,ia,ja,avals,perm,invp,rhs,1,iparm,dparm); fprintf(stdout,"-- FIN INIT PARAMETERS --\n"); init_raff = iparm[IPARM_ITERMAX]; fflush(stdout); /* Passage en mode verbose */ iparm[IPARM_RHS_MAKING] = API_RHS_B; if( !param_int && data_option.empty() ){ iparm[IPARM_MATRIX_VERIFICATION] = API_YES; iparm[IPARM_REFINEMENT] = API_RAF_GMRES; iparm[IPARM_INCOMPLETE] = API_NO; } if( !param_double && data_option.empty()){ dparm[DPARM_EPSILON_REFINEMENT] = 1e-12; dparm[DPARM_EPSILON_MAGN_CTRL] = 1e-32; } SYM = AA.symetrique; cout << "SYM = "<< SYM << endl; // SYMETRIQUE if( SYM == 1 ){ iparm[IPARM_SYM] = API_SYM_YES; //iparm[IPARM_FACTORIZATION] = API_FACT_LDLT; } if( SYM == 0 ){ iparm[IPARM_SYM] = API_SYM_NO; //iparm[IPARM_FACTORIZATION] = API_FACT_LU; } /* Scotch */ fprintf(stdout,"-- Scotch --\n"); fflush(stdout); iparm[IPARM_START_TASK] = API_TASK_ORDERING; iparm[IPARM_END_TASK] = API_TASK_ORDERING; pastix(&pastix_data, MPI_COMM_WORLD, Ncol,ia,ja,avals,perm,invp,rhs,1,iparm,dparm); /* Fax */ fprintf(stdout,"-- Fax --\n"); iparm[IPARM_START_TASK] = API_TASK_SYMBFACT; iparm[IPARM_END_TASK] = API_TASK_SYMBFACT; pastix(&pastix_data, MPI_COMM_WORLD, Ncol,ia,ja,avals,perm,invp,rhs,1,iparm,dparm); /* Blend */ fprintf(stdout,"-- Blend --\n"); iparm[IPARM_START_TASK] = API_TASK_ANALYSE; iparm[IPARM_END_TASK] = API_TASK_ANALYSE; pastix(&pastix_data, MPI_COMM_WORLD, Ncol,ia,ja,avals,perm,invp,rhs,1,iparm,dparm); if( SYM == 1 ){ //iparm[IPARM_SYM] = API_SYM_YES; iparm[IPARM_FACTORIZATION] = API_FACT_LDLT; } if( SYM == 0 ){ //iparm[IPARM_SYM] = API_SYM_NO; iparm[IPARM_FACTORIZATION] = API_FACT_LU; } /* Factorisation */ iparm[IPARM_START_TASK] = API_TASK_NUMFACT; iparm[IPARM_END_TASK] = API_TASK_NUMFACT; gettimeofday(&tv1, NULL); fprintf(stdout,"-- SOPALIN --\n"); pastix(&pastix_data, MPI_COMM_WORLD, Ncol,ia,ja,avals,perm,invp,rhs,1,iparm,dparm); gettimeofday(&tv2, NULL); fprintf(stdout,"Time to call factorization : %ld usec\n", (long)((tv2.tv_sec - tv1.tv_sec ) * 1000000 + tv2.tv_usec - tv1.tv_usec)); for(int ii=0; ii < Ncol+1; ii++) ia[ii] = ia[ii]-1; for(int ii=0; ii < ia[Ncol]-1; ii++) ja[ii] = ja[ii]-1; }