MPI_Comm NM_MPI_com(NumericsMatrix* A) { if (!A || (A && NM_linearSolverParams(A)->mpi_com == MPI_COMM_NULL)) { int myid; int argc = 0; /* C99 requires that argv[argc] == NULL. With openmpi 1.8, we get a * segfault if this is not true */ char *argv0 = NULL; char **argv = &argv0; CHECK_MPI(MPI_Init(&argc, &argv)); CHECK_MPI(MPI_Comm_rank(MPI_COMM_WORLD, &myid)); if (A) { NM_linearSolverParams(A)->mpi_com = MPI_COMM_WORLD; NM_linearSolverParams(A)->mpi_com_init = 1; } } if(A) { return NM_linearSolverParams(A)->mpi_com; } else { return MPI_COMM_WORLD; } }
NM_UMFPACK_WS* NM_UMFPACK_factorize(NumericsMatrix* A) { NumericsSparseLinearSolverParams* params = NM_linearSolverParams(A); if (params->solver_data) { return (NM_UMFPACK_WS*) params->solver_data; } params->solver_data = (NM_UMFPACK_WS*)calloc(1, sizeof(NM_UMFPACK_WS)); NM_UMFPACK_WS* umfpack_ws = (NM_UMFPACK_WS*) params->solver_data; UMFPACK_FN(defaults) (umfpack_ws->control); umfpack_ws->control[UMFPACK_PRL] = verbose; /* TODO UMFPACK_PIVOT_TOLERANCE, UMFPACK_ORDERING, UMFPACK_SCALE * UMFPACK_DROPTOL, UMFPACK_STRATEGY, UMFPACK_IRSTEP*/ CSparseMatrix* C = NM_csc(A); csi status; status = UMFPACK_FN(symbolic) (C->m, C->n, C->p, C->i, C->x, &(umfpack_ws->symbolic), umfpack_ws->control, umfpack_ws->info); if (status) { umfpack_ws->control[UMFPACK_PRL] = 1; UMFPACK_FN(report_status) (umfpack_ws->control, status); return NULL; } status = UMFPACK_FN(numeric) (C->p, C->i, C->x, umfpack_ws->symbolic, &(umfpack_ws->numeric), umfpack_ws->control, umfpack_ws->info); if (status) { umfpack_ws->control[UMFPACK_PRL] = 1; UMFPACK_FN(report_status) (umfpack_ws->control, status); return NULL; } umfpack_ws->wi = (csi*)malloc(C->n * sizeof(csi)); csi size_wd; if (umfpack_ws->control[UMFPACK_IRSTEP] > 0) { size_wd = 5 * C->n; } else { size_wd = C->n; } umfpack_ws->wd = (double*)malloc(size_wd * sizeof(double)); umfpack_ws->x = (double*)malloc(C->n * sizeof(double)); return umfpack_ws; }
void fc3d_nonsmooth_Newton_solvers_solve(fc3d_nonsmooth_Newton_solvers* equation, double* reaction, double* velocity, int* info, SolverOptions* options) { assert(equation); FrictionContactProblem* problem = equation->problem; assert(problem); assert(reaction); assert(velocity); assert(info); assert(options); assert(problem->dimension == 3); assert(options->iparam); assert(options->dparam); assert(problem->q); assert(problem->mu); assert(problem->M); assert(problem->M->matrix0 || problem->M->matrix1 || problem->M->matrix2); assert(!options->iparam[4]); // only host unsigned int problemSize = 3 * problem->numberOfContacts; unsigned int iter = 0; unsigned int itermax = options->iparam[0]; unsigned int erritermax = options->iparam[7]; int nzmax; if (problem->M->storageType == NM_DENSE) { nzmax = problemSize * problemSize; } else { nzmax = options->iparam[3]; } assert(itermax > 0); assert(nzmax > 0); double tolerance = options->dparam[0]; assert(tolerance > 0); if (verbose > 0) printf("------------------------ FC3D - _nonsmooth_Newton_solversSolve - Start with tolerance = %g\n", tolerance); unsigned int _3problemSize = 3 * problemSize; double normq = cblas_dnrm2(problemSize , problem->q , 1); void *buffer; if (!options->dWork) { buffer = malloc((11 * problemSize) * sizeof(double)); // F(1), // tmp1(1), // tmp2(1), // tmp3(1), // A(3), // B(3), rho } else { buffer = options->dWork; } double *F = (double *) buffer; double *tmp1 = (double *) F + problemSize; double *tmp2 = (double *) tmp1 + problemSize; double *tmp3 = (double *) tmp2 + problemSize; double *Ax = tmp3 + problemSize; double *Bx = Ax + _3problemSize; double *rho = Bx + _3problemSize; NumericsMatrix *AWpB, *AWpB_backup; if (!options->dWork) { AWpB = createNumericsMatrix(problem->M->storageType, problem->M->size0, problem->M->size1); AWpB_backup = createNumericsMatrix(problem->M->storageType, problem->M->size0, problem->M->size1); } else { AWpB = (NumericsMatrix*) (rho + problemSize); AWpB_backup = (NumericsMatrix*) (AWpB + sizeof(NumericsMatrix*)); } /* just for allocations */ NM_copy(problem->M, AWpB); if (problem->M->storageType != NM_DENSE) { switch(options->iparam[13]) { case 0: { NM_linearSolverParams(AWpB)->solver = NS_CS_LUSOL; break; } case 1: { NM_linearSolverParams(AWpB)->solver = NS_MUMPS; #ifdef HAVE_MPI assert (options->solverData); if ((MPI_Comm) options->solverData == MPI_COMM_NULL) { options->solverData = NM_MPI_com(MPI_COMM_NULL); } else { NM_MPI_com((MPI_Comm) options->solverData); } #endif break; } default: { numerics_error("fc3d_nonsmooth_Newton_solvers_solve", "Unknown linear solver.\n"); } } } // compute rho here for (unsigned int i = 0; i < problemSize; ++i) rho[i] = options->dparam[3]; // velocity <- M*reaction + qfree cblas_dcopy(problemSize, problem->q, 1, velocity, 1); NM_gemv(1., problem->M, reaction, 1., velocity); double linear_solver_residual=0.0; while (iter++ < itermax) { equation->function(equation->data, problemSize, reaction, velocity, equation->problem->mu, rho, F, Ax, Bx); // AW + B computeAWpB(Ax, problem->M, Bx, AWpB); cblas_dcopy_msan(problemSize, F, 1, tmp1, 1); cblas_dscal(problemSize, -1., tmp1, 1); /* Solve: AWpB X = -F */ NM_copy(AWpB, AWpB_backup); int lsi = NM_gesv(AWpB, tmp1); /* NM_copy needed here */ NM_copy(AWpB_backup, AWpB); if (lsi) { if (verbose > 0) { numerics_warning("fc3d_nonsmooth_Newton_solvers_solve", "warning! linear solver exit with code = %d\n", lsi); } } if (verbose > 0) { cblas_dcopy_msan(problemSize, F, 1, tmp3, 1); NM_gemv(1., AWpB, tmp1, 1., tmp3); linear_solver_residual = cblas_dnrm2(problemSize, tmp3, 1); /* fprintf(stderr, "fc3d esolve: linear equation residual = %g\n", */ /* cblas_dnrm2(problemSize, tmp3, 1)); */ /* for the component wise scaled residual: cf mumps & * http://www.netlib.org/lapack/lug/node81.html */ } // line search double alpha = 1; int info_ls = 0; cblas_dcopy_msan(problemSize, tmp1, 1, tmp3, 1); switch (options->iparam[11]) { case -1: /* without line search */ info_ls = 1; break; case 0: /* Goldstein Price */ info_ls = globalLineSearchGP(equation, reaction, velocity, problem->mu, rho, F, Ax, Bx, problem->M, problem->q, AWpB, tmp1, tmp2, &alpha, options->iparam[12]); break; case 1: /* FBLSA */ info_ls = frictionContactFBLSA(equation, reaction, velocity, problem->mu, rho, F, Ax, Bx, problem->M, problem->q, AWpB, tmp1, tmp2, &alpha, options->iparam[12]); break; default: { numerics_error("fc3d_nonsmooth_Newton_solvers_solve", "Unknown line search option.\n"); } } if (!info_ls) // tmp2 should contains the reaction iterate of the line search // for GP this should be the same as cblas_daxpy(problemSize, alpha, tmp1, 1, reaction, 1); cblas_dcopy(problemSize, tmp2, 1, reaction, 1); else cblas_daxpy(problemSize, 1., tmp3, 1., reaction, 1); // velocity <- M*reaction + qfree cblas_dcopy(problemSize, problem->q, 1, velocity, 1); NM_gemv(1., problem->M, reaction, 1., velocity); options->dparam[1] = INFINITY; if (!(iter % erritermax)) { fc3d_compute_error(problem, reaction, velocity, // fc3d_FischerBurmeister_compute_error(problem, reaction, velocity, tolerance, options, normq, &(options->dparam[1])); DEBUG_EXPR_WE(equation->function(equation->data, problemSize, reaction, velocity, equation->problem->mu, rho, F, NULL, NULL)); DEBUG_EXPR_WE(assert((cblas_dnrm2(problemSize, F, 1) / (1 + cblas_dnrm2(problemSize, problem->q, 1))) <= (10 * options->dparam[1] + 1e-15))); } if (verbose > 0) { equation->function(equation->data, problemSize, reaction, velocity, equation->problem->mu, rho, F, NULL, NULL); printf(" ---- fc3d_nonsmooth_Newton_solvers_solve: iteration %d : , linear solver residual =%g, residual=%g, ||F||=%g\n", iter, linear_solver_residual, options->dparam[1],cblas_dnrm2(problemSize, F, 1)); } if (options->callback) { options->callback->collectStatsIteration(options->callback->env, problemSize, reaction, velocity, options->dparam[1], NULL); } if (isnan(options->dparam[1])) { if (verbose > 0) { printf(" fc3d_nonsmooth_Newton_solvers_solve: iteration %d : computed residual is not a number, stop.\n", iter); } info[0] = 2; break; } if (options->dparam[1] < tolerance) { info[0] = 0; break; } } if (verbose > 0) { if (!info[0]) printf("------------------------ FC3D - NSN - convergence after %d iterations, residual : %g < %g \n", iter, options->dparam[1],tolerance); else { printf("------------------------ FC3D - NSN - no convergence after %d iterations, residual : %g < %g \n", iter, options->dparam[1], tolerance); } } options->iparam[SICONOS_IPARAM_ITER_DONE] = iter; if (!options->dWork) { assert(buffer); free(buffer); options->dWork = NULL; } else { assert(buffer == options->dWork); } if (!options->dWork) { freeNumericsMatrix(AWpB); freeNumericsMatrix(AWpB_backup); free(AWpB); free(AWpB_backup); } if (verbose > 0) printf("------------------------ FC3D - NSN - End\n"); }
int NM_gesv(NumericsMatrix* A, double *b) { assert(A->size0 == A->size1); int info = 1; switch (A->storageType) { case NM_DENSE: { assert(A->matrix0); DGESV(A->size0, 1, A->matrix0, A->size0, NM_iWork(A, A->size0), b, A->size0, &info); break; } case NM_SPARSE_BLOCK: /* sparse block -> triplet -> csc */ case NM_SPARSE: { switch (NM_linearSolverParams(A)->solver) { case NS_CS_LUSOL: info = !cs_lusol(1, NM_csc(A), b, DBL_EPSILON); break; #ifdef WITH_MUMPS case NS_MUMPS: { /* the mumps instance is initialized (call with job=-1) */ DMUMPS_STRUC_C* mumps_id = NM_MUMPS_id(A); mumps_id->rhs = b; mumps_id->job = 6; /* compute the solution */ dmumps_c(mumps_id); /* clean the mumps instance */ mumps_id->job = -2; dmumps_c(mumps_id); info = mumps_id->info[0]; if (info > 0) { if (verbose > 0) { printf("NM_gesv: MUMPS fails : info(1)=%d, info(2)=%d\n", info, mumps_id->info[1]); } } if (verbose > 1) { printf("MUMPS : condition number %g\n", mumps_id->rinfog[9]); printf("MUMPS : component wise scaled residual %g\n", mumps_id->rinfog[6]); printf("MUMPS : \n"); } /* Here we free mumps_id ... */ free(NM_linearSolverParams(A)->solver_data); NM_linearSolverParams(A)->solver_data = NULL; break; } #endif default: { fprintf(stderr, "NM_gesv: unknown sparse linearsolver : %d\n", NM_linearSolverParams(A)->solver); exit(EXIT_FAILURE); } } break; } default: assert (0 && "NM_gesv unknown storageType"); } /* some time we cannot find a solution to a linear system, and its fine, for * instance with the minFBLSA. Therefore, we should not check here for * problems, but the calling function has to check the return code.*/ // CHECK_RETURN(info); return info; }
DMUMPS_STRUC_C* NM_MUMPS_id(NumericsMatrix* A) { NumericsSparseLinearSolverParams* params = NM_linearSolverParams(A); if (!params->solver_data) { params->solver_data = malloc(sizeof(DMUMPS_STRUC_C)); DMUMPS_STRUC_C* mumps_id = (DMUMPS_STRUC_C*) params->solver_data; // Initialize a MUMPS instance. Use MPI_COMM_WORLD. mumps_id->job = JOB_INIT; mumps_id->par = 1; mumps_id->sym = 0; if (NM_MPI_com(A) == MPI_COMM_WORLD) { mumps_id->comm_fortran = USE_COMM_WORLD; } else { mumps_id->comm_fortran = MPI_Comm_c2f(NM_MPI_com(A)); } dmumps_c(mumps_id); if (verbose == 1) { mumps_id->ICNTL(1) = -1; // Error messages, standard output stream. mumps_id->ICNTL(2) = -1; // Diagnostics, standard output stream. mumps_id->ICNTL(3) = -1; // Global infos, standard output stream. mumps_id->ICNTL(11) = 1; // Error analysis } else if (verbose == 2) { mumps_id->ICNTL(1) = -1; // Error messages, standard output stream. mumps_id->ICNTL(2) = -1; // Diagnostics, standard output stream. mumps_id->ICNTL(3) = 6; // Global infos, standard output stream. // mumps_id->ICNTL(4) = 4; // Errors, warnings and information on // input, output parameters printed. // mumps_id->ICNTL(10) = 1; // One step of iterative refinment mumps_id->ICNTL(11) = 1; // Error analysis } else if (verbose >= 3) { mumps_id->ICNTL(1) = 6; // Error messages, standard output stream. mumps_id->ICNTL(2) = 6; // Diagnostics, standard output stream. mumps_id->ICNTL(3) = 6; // Global infos, standard output stream. // mumps_id->ICNTL(4) = 4; // Errors, warnings and information on // input, output parameters printed. // mumps_id->ICNTL(10) = 1; // One step of iterative refinment mumps_id->ICNTL(11) = 1; // Error analysis } else { mumps_id->ICNTL(1) = -1; mumps_id->ICNTL(2) = -1; mumps_id->ICNTL(3) = -1; } mumps_id->ICNTL(24) = 1; // Null pivot row detection see also CNTL(3) & CNTL(5) // ok for a cube on a plane & four contact points // computeAlartCurnierSTD != generated in this case... //mumps_id->CNTL(3) = ...; //mumps_id->CNTL(5) = ...; } DMUMPS_STRUC_C* mumps_id = (DMUMPS_STRUC_C*) params->solver_data; mumps_id->n = (int) NM_triplet(A)->n; mumps_id->irn = NM_MUMPS_irn(A); mumps_id->jcn = NM_MUMPS_jcn(A); int nz; if (NM_sparse(A)->triplet) { nz = (int) NM_sparse(A)->triplet->nz; mumps_id->nz = nz; mumps_id->a = NM_sparse(A)->triplet->x; } else { nz = NM_linearSolverParams(A)->iWork[2 * NM_csc(A)->nzmax]; mumps_id->nz = nz; mumps_id->a = NM_sparse(A)->csc->x; } return (DMUMPS_STRUC_C*) params->solver_data; }