void ShapeImprovementWrapper::run_wrapper( Mesh* mesh, ParallelMesh* pmesh, MeshDomain* domain, Settings* settings, QualityAssessor* qa, MsqError& err ) { // Define an untangler UntangleBetaQualityMetric untangle_metric( untBeta ); LPtoPTemplate untangle_func( 2, &untangle_metric ); ConjugateGradient untangle_global( &untangle_func ); TerminationCriterion untangle_inner, untangle_outer; untangle_global.use_global_patch(); untangle_inner.add_absolute_quality_improvement( 0.0 ); untangle_inner.add_absolute_successive_improvement( successiveEps ); untangle_outer.add_iteration_limit( 1 ); untangle_global.set_inner_termination_criterion( &untangle_inner ); untangle_global.set_outer_termination_criterion( &untangle_outer ); // define shape improver IdealWeightInverseMeanRatio inverse_mean_ratio; inverse_mean_ratio.set_averaging_method( QualityMetric::LINEAR ); LPtoPTemplate obj_func( 2, &inverse_mean_ratio ); FeasibleNewton feas_newt( &obj_func ); TerminationCriterion term_inner, term_outer; feas_newt.use_global_patch(); qa->add_quality_assessment( &inverse_mean_ratio ); term_inner.add_absolute_gradient_L2_norm( gradNorm ); term_inner.add_relative_successive_improvement( successiveEps ); term_outer.add_iteration_limit( pmesh ? parallelIterations : 1 ); feas_newt.set_inner_termination_criterion( &term_inner ); feas_newt.set_outer_termination_criterion( &term_outer ); // Apply CPU time limit to untangler if (maxTime > 0.0) untangle_inner.add_cpu_time( maxTime ); // Run untangler InstructionQueue q1; Timer totalTimer; q1.set_master_quality_improver( &untangle_global, err ); MSQ_ERRRTN(err); q1.add_quality_assessor( qa, err ); MSQ_ERRRTN(err); q1.run_common( mesh, pmesh, domain, settings, err ); MSQ_ERRRTN(err); // If limited by CPU time, limit next step to remaning time if (maxTime > 0.0) { double remaining = maxTime - totalTimer.since_birth(); if (remaining <= 0.0 ){ MSQ_DBGOUT(2) << "Optimization is terminating without perfoming shape improvement." << std::endl; remaining = 0.0; } term_inner.add_cpu_time( remaining ); } // Run shape improver InstructionQueue q2; q2.set_master_quality_improver( &feas_newt, err ); MSQ_ERRRTN(err); q2.add_quality_assessor( qa, err ); MSQ_ERRRTN(err); q2.run_common( mesh, pmesh, domain, settings, err ); MSQ_ERRRTN(err); }
void PMMParallelReferenceMeshSmoother::run_wrapper( Mesh* mesh, ParallelMesh* pmesh, MeshDomain* domain, Settings* settings, QualityAssessor* qa, MsqError& err ) { std::cout << "\nP[" << Mesquite::get_parallel_rank() << "] tmp srk PMMParallelReferenceMeshSmoother innerIter= " << innerIter << " parallelIterations= " << parallelIterations << std::endl; //if (!get_parallel_rank()) std::cout << "\nP[" << get_parallel_rank() << "] tmp srk PMMParallelReferenceMeshSmoother: running shape improver... \n" << std::endl; PerceptMesquiteMesh *pmm = dynamic_cast<PerceptMesquiteMesh *>(mesh); PerceptMesh *eMesh = pmm->getPerceptMesh(); m_pmm= pmm; m_eMesh = eMesh; m_num_nodes = m_eMesh->get_number_nodes(); print_comm_list(*eMesh->get_bulk_data(), false); stk::mesh::FieldBase *coord_field = eMesh->get_coordinates_field(); stk::mesh::FieldBase *coord_field_current = coord_field; stk::mesh::FieldBase *coord_field_projected = eMesh->get_field("coordinates_N"); stk::mesh::FieldBase *coord_field_original = eMesh->get_field("coordinates_NM1"); stk::mesh::FieldBase *coord_field_lagged = eMesh->get_field("coordinates_lagged"); m_coord_field_original = coord_field_original; m_coord_field_projected = coord_field_projected; m_coord_field_lagged = coord_field_lagged; m_coord_field_current = coord_field_current; eMesh->copy_field(coord_field_lagged, coord_field_original); // untangle PMMSmootherMetricUntangle untangle_metric(eMesh); // shape-size-orient smooth PMMSmootherMetricShapeSizeOrient shape_metric(eMesh); // shape PMMSmootherMetricShapeB1 shape_b1_metric(eMesh); // scaled jacobian PMMSmootherMetricScaledJacobianElemental scaled_jac_metric(eMesh); // scaled jacobian - nodal PMMSmootherMetricScaledJacobianNodal scaled_jac_metric_nodal(eMesh); //double omegas[] = {0.0, 0.001, 0.01, 0.1, 0.2, 0.4, 0.6, 0.8, 1.0}; //double omegas[] = {0.001, 1.0}; //double omegas[] = { 0.001, 0.01, 0.02, 0.03, 0.04, 0.05, 0.06, 0.07, 0.08, 0.09, 0.1, 0.2, 0.4, 0.45,0.46,0.47,0.48,0.49,0.5,0.52,0.54,0.56,0.59, 0.6, 0.8, 1.0}; double omegas[] = { 1.0}; //double omegas[] = {0.0, 0.001, 0.01, 0.1, 0.11, 0.12, 0.13, 0.14, 0.15, 0.16, 0.17, 0.18, 0.19, 0.2, 0.4, 0.6, 0.8, 1.0}; int nomega = sizeof(omegas)/sizeof(omegas[0]); for (int outer = 0; outer < nomega; outer++) { double omega = (outer < nomega ? omegas[outer] : 1.0); m_omega = omega; m_omega_prev = omega; if (outer > 0) m_omega_prev = omegas[outer-1]; // set current state and evaluate mesh validity (current = omega*project + (1-omega)*original) eMesh->nodal_field_axpbypgz(omega, coord_field_projected, (1.0-omega), coord_field_original, 0.0, coord_field_current); int num_invalid = PMMParallelShapeImprover::parallel_count_invalid_elements(m_eMesh); if (!get_parallel_rank()) std::cout << "\ntmp srk PMMParallelReferenceMeshSmoother num_invalid current= " << num_invalid << " for outer_iter= " << outer << " omega= " << omega << (num_invalid ? " WARNING: invalid elements exist before Mesquite smoothing" : " OK") << std::endl; //if (num_invalid) return; m_num_invalid = num_invalid; m_untangled = (m_num_invalid == 0); int iter_all=0; int do_anim = 0; // = frequency of anim writes if (do_anim) { eMesh->save_as("anim_all."+toString(iter_all)+".e"); } for (int stage = 0; stage < 2; stage++) { m_stage = stage; if (stage==0) { m_metric = &untangle_metric; //m_metric = &scaled_jac_metric_nodal; } else { int num_invalid_1 = PMMParallelShapeImprover::parallel_count_invalid_elements(m_eMesh); VERIFY_OP_ON(num_invalid_1, ==, 0, "Invalid elements exist for start of stage 2, aborting"); //m_metric = &shape_metric; m_metric = &shape_b1_metric; //m_metric = &scaled_jac_metric; } for (int iter = 0; iter < innerIter; ++iter, ++iter_all) { m_iter = iter; int num_invalid_0 = PMMParallelShapeImprover::parallel_count_invalid_elements(m_eMesh); m_num_invalid = num_invalid_0; // if (!get_parallel_rank() && num_invalid_0) // std::cout << "\ntmp srk PMMParallelReferenceMeshSmoother num_invalid current= " << num_invalid_0 // << (num_invalid ? " WARNING: invalid elements exist before Mesquite smoothing" : "OK") // << std::endl; m_global_metric = run_one_iteration(mesh, domain, err); sync_fields(iter); num_invalid_0 = PMMParallelShapeImprover::parallel_count_invalid_elements(m_eMesh); m_num_invalid = num_invalid_0; bool conv = check_convergence(); if (!get_parallel_rank()) { std::cout << "P[" << get_parallel_rank() << "] " << "tmp srk iter= " << iter << " dmax= " << m_dmax << " m_dnew= " << m_dnew << " m_d0= " << m_d0 << " m_alpha= " << m_alpha << " m_grad_norm= " << m_grad_norm << " m_scaled_grad_norm = " << m_scaled_grad_norm << " num_invalid= " << num_invalid_0 << " m_global_metric= " << m_global_metric << " m_untangled= " << m_untangled << std::endl; } if (do_anim) { eMesh->save_as("iter_"+toString(outer)+"_"+toString(stage)+"."+toString(iter+1)+".e"); if (iter_all % do_anim == 0) eMesh->save_as("anim_all."+toString(iter_all+1)+".e"); } if (!m_untangled && m_num_invalid == 0) { m_untangled = true; } if (conv && m_untangled) break; //if (iter == 5) break; //if (iter == 0) exit(1); } eMesh->save_as("outer_iter_"+toString(outer)+"_"+toString(stage)+"_mesh.e"); } eMesh->copy_field(coord_field_lagged, coord_field); } //if (!get_parallel_rank()) MPI_Barrier( MPI_COMM_WORLD ); std::cout << "\nP[" << get_parallel_rank() << "] tmp srk PMMParallelReferenceMeshSmoother: running shape improver... done \n" << std::endl; MSQ_ERRRTN(err); }