1674ae819SStefano Zampini #include "bddc.h" 2674ae819SStefano Zampini #include "bddcprivate.h" 3674ae819SStefano Zampini #include <petscblaslapack.h> 4674ae819SStefano Zampini 5674ae819SStefano Zampini #undef __FUNCT__ 6674ae819SStefano Zampini #define __FUNCT__ "PCBDDCResetCustomization" 7674ae819SStefano Zampini PetscErrorCode PCBDDCResetCustomization(PC pc) 8674ae819SStefano Zampini { 9674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 10674ae819SStefano Zampini PetscInt i; 11674ae819SStefano Zampini PetscErrorCode ierr; 12674ae819SStefano Zampini 13674ae819SStefano Zampini PetscFunctionBegin; 14674ae819SStefano Zampini ierr = PCBDDCGraphResetCSR(pcbddc->mat_graph);CHKERRQ(ierr); 15674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr); 16674ae819SStefano Zampini ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr); 17674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 18674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 19674ae819SStefano Zampini for (i=0;i<pcbddc->n_ISForDofs;i++) { 20674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 21674ae819SStefano Zampini } 22674ae819SStefano Zampini ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 23674ae819SStefano Zampini PetscFunctionReturn(0); 24674ae819SStefano Zampini } 25674ae819SStefano Zampini 26674ae819SStefano Zampini #undef __FUNCT__ 27674ae819SStefano Zampini #define __FUNCT__ "PCBDDCResetTopography" 28674ae819SStefano Zampini PetscErrorCode PCBDDCResetTopography(PC pc) 29674ae819SStefano Zampini { 30674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 31674ae819SStefano Zampini PetscErrorCode ierr; 32674ae819SStefano Zampini 33674ae819SStefano Zampini PetscFunctionBegin; 34674ae819SStefano Zampini ierr = MatDestroy(&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr); 35674ae819SStefano Zampini ierr = MatDestroy(&pcbddc->ConstraintMatrix);CHKERRQ(ierr); 36674ae819SStefano Zampini ierr = PCBDDCGraphReset(pcbddc->mat_graph);CHKERRQ(ierr); 37674ae819SStefano Zampini PetscFunctionReturn(0); 38674ae819SStefano Zampini } 39674ae819SStefano Zampini 40674ae819SStefano Zampini #undef __FUNCT__ 41674ae819SStefano Zampini #define __FUNCT__ "PCBDDCResetSolvers" 42674ae819SStefano Zampini PetscErrorCode PCBDDCResetSolvers(PC pc) 43674ae819SStefano Zampini { 44674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 45674ae819SStefano Zampini PetscErrorCode ierr; 46674ae819SStefano Zampini 47674ae819SStefano Zampini PetscFunctionBegin; 48674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 49674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 50674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->coarse_vec);CHKERRQ(ierr); 51674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->coarse_rhs);CHKERRQ(ierr); 52674ae819SStefano Zampini ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr); 53674ae819SStefano Zampini ierr = MatDestroy(&pcbddc->coarse_mat);CHKERRQ(ierr); 54674ae819SStefano Zampini ierr = MatDestroy(&pcbddc->coarse_phi_B);CHKERRQ(ierr); 55674ae819SStefano Zampini ierr = MatDestroy(&pcbddc->coarse_phi_D);CHKERRQ(ierr); 5615aaf578SStefano Zampini ierr = MatDestroy(&pcbddc->coarse_psi_B);CHKERRQ(ierr); 5715aaf578SStefano Zampini ierr = MatDestroy(&pcbddc->coarse_psi_D);CHKERRQ(ierr); 58674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->vec1_P);CHKERRQ(ierr); 59674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->vec1_C);CHKERRQ(ierr); 60674ae819SStefano Zampini ierr = MatDestroy(&pcbddc->local_auxmat1);CHKERRQ(ierr); 61674ae819SStefano Zampini ierr = MatDestroy(&pcbddc->local_auxmat2);CHKERRQ(ierr); 62674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->vec1_R);CHKERRQ(ierr); 63674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->vec2_R);CHKERRQ(ierr); 64674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->vec4_D);CHKERRQ(ierr); 65674ae819SStefano Zampini ierr = VecScatterDestroy(&pcbddc->R_to_B);CHKERRQ(ierr); 66674ae819SStefano Zampini ierr = VecScatterDestroy(&pcbddc->R_to_D);CHKERRQ(ierr); 67674ae819SStefano Zampini ierr = VecScatterDestroy(&pcbddc->coarse_loc_to_glob);CHKERRQ(ierr); 68674ae819SStefano Zampini ierr = PetscFree(pcbddc->local_primal_indices);CHKERRQ(ierr); 69674ae819SStefano Zampini ierr = PetscFree(pcbddc->replicated_local_primal_indices);CHKERRQ(ierr); 70674ae819SStefano Zampini ierr = PetscFree(pcbddc->replicated_local_primal_values);CHKERRQ(ierr); 71674ae819SStefano Zampini ierr = PetscFree(pcbddc->local_primal_displacements);CHKERRQ(ierr); 72674ae819SStefano Zampini ierr = PetscFree(pcbddc->local_primal_sizes);CHKERRQ(ierr); 73674ae819SStefano Zampini PetscFunctionReturn(0); 74674ae819SStefano Zampini } 75674ae819SStefano Zampini 76674ae819SStefano Zampini #undef __FUNCT__ 77aa0d41d4SStefano Zampini #define __FUNCT__ "PCBDDCSetUpLocalMatrices" 78aa0d41d4SStefano Zampini PetscErrorCode PCBDDCSetUpLocalMatrices(PC pc) 79aa0d41d4SStefano Zampini { 80aa0d41d4SStefano Zampini PC_IS* pcis = (PC_IS*)(pc->data); 81aa0d41d4SStefano Zampini PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 82aa0d41d4SStefano Zampini Mat_IS* matis = (Mat_IS*)pc->pmat->data; 83aa0d41d4SStefano Zampini /* manage repeated solves */ 84aa0d41d4SStefano Zampini MatReuse reuse; 85aa0d41d4SStefano Zampini MatStructure matstruct; 86aa0d41d4SStefano Zampini PetscErrorCode ierr; 87aa0d41d4SStefano Zampini 88aa0d41d4SStefano Zampini PetscFunctionBegin; 89aa0d41d4SStefano Zampini /* get mat flags */ 90aa0d41d4SStefano Zampini ierr = PCGetOperators(pc,NULL,NULL,&matstruct);CHKERRQ(ierr); 91aa0d41d4SStefano Zampini reuse = MAT_INITIAL_MATRIX; 92aa0d41d4SStefano Zampini if (pc->setupcalled) { 93aa0d41d4SStefano Zampini /* when matstruct is SAME_PRECONDITIONER, we shouldn't be here */ 94aa0d41d4SStefano Zampini if (matstruct == SAME_PRECONDITIONER) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"This should not happen"); 95aa0d41d4SStefano Zampini if (matstruct == SAME_NONZERO_PATTERN) { 96aa0d41d4SStefano Zampini reuse = MAT_REUSE_MATRIX; 97aa0d41d4SStefano Zampini } else { 98aa0d41d4SStefano Zampini reuse = MAT_INITIAL_MATRIX; 99aa0d41d4SStefano Zampini } 100aa0d41d4SStefano Zampini } 101aa0d41d4SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 102aa0d41d4SStefano Zampini ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr); 103aa0d41d4SStefano Zampini ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr); 104aa0d41d4SStefano Zampini ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr); 105aa0d41d4SStefano Zampini ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr); 106aa0d41d4SStefano Zampini ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 107aa0d41d4SStefano Zampini } 108aa0d41d4SStefano Zampini 109aa0d41d4SStefano Zampini /* transform local matrices if needed */ 110aa0d41d4SStefano Zampini if (pcbddc->use_change_of_basis) { 111aa0d41d4SStefano Zampini Mat change_mat_all; 112aa0d41d4SStefano Zampini PetscScalar *row_cmat_values; 113aa0d41d4SStefano Zampini PetscInt *row_cmat_indices; 114aa0d41d4SStefano Zampini PetscInt *nnz,*is_indices,*temp_indices; 115aa0d41d4SStefano Zampini PetscInt i,j,k,n_D,n_B; 116aa0d41d4SStefano Zampini 117aa0d41d4SStefano Zampini /* Get Non-overlapping dimensions */ 118aa0d41d4SStefano Zampini n_B = pcis->n_B; 119aa0d41d4SStefano Zampini n_D = pcis->n-n_B; 120aa0d41d4SStefano Zampini 121aa0d41d4SStefano Zampini /* compute nonzero structure of change of basis on all local nodes */ 122aa0d41d4SStefano Zampini ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 123aa0d41d4SStefano Zampini ierr = ISGetIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 124aa0d41d4SStefano Zampini for (i=0;i<n_D;i++) nnz[is_indices[i]] = 1; 125aa0d41d4SStefano Zampini ierr = ISRestoreIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 126aa0d41d4SStefano Zampini ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 127aa0d41d4SStefano Zampini k=1; 128aa0d41d4SStefano Zampini for (i=0;i<n_B;i++) { 129aa0d41d4SStefano Zampini ierr = MatGetRow(pcbddc->ChangeOfBasisMatrix,i,&j,NULL,NULL);CHKERRQ(ierr); 130aa0d41d4SStefano Zampini nnz[is_indices[i]]=j; 131aa0d41d4SStefano Zampini if (k < j) k = j; 132aa0d41d4SStefano Zampini ierr = MatRestoreRow(pcbddc->ChangeOfBasisMatrix,i,&j,NULL,NULL);CHKERRQ(ierr); 133aa0d41d4SStefano Zampini } 134aa0d41d4SStefano Zampini ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 135aa0d41d4SStefano Zampini /* assemble change of basis matrix on the whole set of local dofs */ 136aa0d41d4SStefano Zampini ierr = PetscMalloc(k*sizeof(PetscInt),&temp_indices);CHKERRQ(ierr); 137aa0d41d4SStefano Zampini ierr = MatCreate(PETSC_COMM_SELF,&change_mat_all);CHKERRQ(ierr); 138aa0d41d4SStefano Zampini ierr = MatSetSizes(change_mat_all,pcis->n,pcis->n,pcis->n,pcis->n);CHKERRQ(ierr); 139aa0d41d4SStefano Zampini ierr = MatSetType(change_mat_all,MATSEQAIJ);CHKERRQ(ierr); 140aa0d41d4SStefano Zampini ierr = MatSeqAIJSetPreallocation(change_mat_all,0,nnz);CHKERRQ(ierr); 141aa0d41d4SStefano Zampini ierr = ISGetIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 142aa0d41d4SStefano Zampini for (i=0;i<n_D;i++) { 143aa0d41d4SStefano Zampini ierr = MatSetValue(change_mat_all,is_indices[i],is_indices[i],1.0,INSERT_VALUES);CHKERRQ(ierr); 144aa0d41d4SStefano Zampini } 145aa0d41d4SStefano Zampini ierr = ISRestoreIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 146aa0d41d4SStefano Zampini ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 147aa0d41d4SStefano Zampini for (i=0;i<n_B;i++) { 148aa0d41d4SStefano Zampini ierr = MatGetRow(pcbddc->ChangeOfBasisMatrix,i,&j,(const PetscInt**)&row_cmat_indices,(const PetscScalar**)&row_cmat_values);CHKERRQ(ierr); 149aa0d41d4SStefano Zampini for (k=0; k<j; k++) temp_indices[k]=is_indices[row_cmat_indices[k]]; 150aa0d41d4SStefano Zampini ierr = MatSetValues(change_mat_all,1,&is_indices[i],j,temp_indices,row_cmat_values,INSERT_VALUES);CHKERRQ(ierr); 151aa0d41d4SStefano Zampini ierr = MatRestoreRow(pcbddc->ChangeOfBasisMatrix,i,&j,(const PetscInt**)&row_cmat_indices,(const PetscScalar**)&row_cmat_values);CHKERRQ(ierr); 152aa0d41d4SStefano Zampini } 153aa0d41d4SStefano Zampini ierr = MatAssemblyBegin(change_mat_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 154aa0d41d4SStefano Zampini ierr = MatAssemblyEnd(change_mat_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 155aa0d41d4SStefano Zampini /* TODO: HOW TO WORK WITH BAIJ? PtAP not provided */ 156aa0d41d4SStefano Zampini ierr = MatGetBlockSize(matis->A,&i);CHKERRQ(ierr); 157aa0d41d4SStefano Zampini if (i==1) { 158aa0d41d4SStefano Zampini ierr = MatPtAP(matis->A,change_mat_all,reuse,2.0,&pcbddc->local_mat);CHKERRQ(ierr); 159aa0d41d4SStefano Zampini } else { 160aa0d41d4SStefano Zampini Mat work_mat; 161aa0d41d4SStefano Zampini ierr = MatConvert(matis->A,MATSEQAIJ,MAT_INITIAL_MATRIX,&work_mat);CHKERRQ(ierr); 162aa0d41d4SStefano Zampini ierr = MatPtAP(work_mat,change_mat_all,reuse,2.0,&pcbddc->local_mat);CHKERRQ(ierr); 163aa0d41d4SStefano Zampini ierr = MatDestroy(&work_mat);CHKERRQ(ierr); 164aa0d41d4SStefano Zampini } 165aa0d41d4SStefano Zampini ierr = MatDestroy(&change_mat_all);CHKERRQ(ierr); 166aa0d41d4SStefano Zampini ierr = PetscFree(nnz);CHKERRQ(ierr); 167aa0d41d4SStefano Zampini ierr = PetscFree(temp_indices);CHKERRQ(ierr); 168aa0d41d4SStefano Zampini } else { 169aa0d41d4SStefano Zampini /* without change of basis, the local matrix is unchanged */ 170aa0d41d4SStefano Zampini if (!pcbddc->local_mat) { 171aa0d41d4SStefano Zampini ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 172aa0d41d4SStefano Zampini pcbddc->local_mat = matis->A; 173aa0d41d4SStefano Zampini } 174aa0d41d4SStefano Zampini } 175aa0d41d4SStefano Zampini 176aa0d41d4SStefano Zampini /* get submatrices */ 177aa0d41d4SStefano Zampini ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_I_local,reuse,&pcis->A_II);CHKERRQ(ierr); 178aa0d41d4SStefano Zampini ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr); 179aa0d41d4SStefano Zampini ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr); 180aa0d41d4SStefano Zampini ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,reuse,&pcis->A_BB);CHKERRQ(ierr); 181aa0d41d4SStefano Zampini PetscFunctionReturn(0); 182aa0d41d4SStefano Zampini } 183aa0d41d4SStefano Zampini 184aa0d41d4SStefano Zampini #undef __FUNCT__ 185*a64d13efSStefano Zampini #define __FUNCT__ "PCBDDCSetUpLocalScatters" 186*a64d13efSStefano Zampini PetscErrorCode PCBDDCSetUpLocalScatters(PC pc,IS* is_R_local_n) 187*a64d13efSStefano Zampini { 188*a64d13efSStefano Zampini PC_IS* pcis = (PC_IS*)(pc->data); 189*a64d13efSStefano Zampini PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 190*a64d13efSStefano Zampini IS is_R_local,is_aux1,is_aux2; 191*a64d13efSStefano Zampini PetscInt *vertices,*aux_array1,*aux_array2,*is_indices,*idx_R_local; 192*a64d13efSStefano Zampini PetscInt n_vertices,n_constraints,i,j,n_R,n_D,n_B; 193*a64d13efSStefano Zampini PetscBool *array_bool; 194*a64d13efSStefano Zampini PetscErrorCode ierr; 195*a64d13efSStefano Zampini 196*a64d13efSStefano Zampini PetscFunctionBegin; 197*a64d13efSStefano Zampini /* Set Non-overlapping dimensions */ 198*a64d13efSStefano Zampini n_B = pcis->n_B; n_D = pcis->n - n_B; 199*a64d13efSStefano Zampini /* get vertex indices from constraint matrix */ 200*a64d13efSStefano Zampini ierr = PCBDDCGetPrimalVerticesLocalIdx(pc,&n_vertices,&vertices);CHKERRQ(ierr); 201*a64d13efSStefano Zampini /* Set number of constraints */ 202*a64d13efSStefano Zampini n_constraints = pcbddc->local_primal_size-n_vertices; 203*a64d13efSStefano Zampini /* Dohrmann's notation: dofs splitted in R (Remaining: all dofs but the vertices) and V (Vertices) */ 204*a64d13efSStefano Zampini ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&array_bool);CHKERRQ(ierr); 205*a64d13efSStefano Zampini for (i=0;i<pcis->n;i++) array_bool[i] = PETSC_TRUE; 206*a64d13efSStefano Zampini for (i=0;i<n_vertices;i++) array_bool[vertices[i]] = PETSC_FALSE; 207*a64d13efSStefano Zampini ierr = PetscMalloc((pcis->n-n_vertices)*sizeof(PetscInt),&idx_R_local);CHKERRQ(ierr); 208*a64d13efSStefano Zampini for (i=0, n_R=0; i<pcis->n; i++) { 209*a64d13efSStefano Zampini if (array_bool[i]) { 210*a64d13efSStefano Zampini idx_R_local[n_R] = i; 211*a64d13efSStefano Zampini n_R++; 212*a64d13efSStefano Zampini } 213*a64d13efSStefano Zampini } 214*a64d13efSStefano Zampini ierr = PetscFree(vertices);CHKERRQ(ierr); 215*a64d13efSStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,n_R,idx_R_local,PETSC_OWN_POINTER,&is_R_local);CHKERRQ(ierr); 216*a64d13efSStefano Zampini 217*a64d13efSStefano Zampini /* print some info if requested */ 218*a64d13efSStefano Zampini if (pcbddc->dbg_flag) { 219*a64d13efSStefano Zampini ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"--------------------------------------------------\n");CHKERRQ(ierr); 220*a64d13efSStefano Zampini ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 221*a64d13efSStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d local dimensions\n",PetscGlobalRank);CHKERRQ(ierr); 222*a64d13efSStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"local_size = %d, dirichlet_size = %d, boundary_size = %d\n",pcis->n,n_D,n_B);CHKERRQ(ierr); 223*a64d13efSStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"r_size = %d, v_size = %d, constraints = %d, local_primal_size = %d\n",n_R,n_vertices,n_constraints,pcbddc->local_primal_size);CHKERRQ(ierr); 224*a64d13efSStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"pcbddc->n_vertices = %d, pcbddc->n_constraints = %d\n",pcbddc->n_vertices,pcbddc->n_constraints);CHKERRQ(ierr); 225*a64d13efSStefano Zampini ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 226*a64d13efSStefano Zampini } 227*a64d13efSStefano Zampini 228*a64d13efSStefano Zampini /* VecScatters pcbddc->R_to_B and (optionally) pcbddc->R_to_D */ 229*a64d13efSStefano Zampini ierr = PetscMalloc((pcis->n_B-n_vertices)*sizeof(PetscInt),&aux_array1);CHKERRQ(ierr); 230*a64d13efSStefano Zampini ierr = PetscMalloc((pcis->n_B-n_vertices)*sizeof(PetscInt),&aux_array2);CHKERRQ(ierr); 231*a64d13efSStefano Zampini ierr = ISGetIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 232*a64d13efSStefano Zampini for (i=0; i<n_D; i++) array_bool[is_indices[i]] = PETSC_FALSE; 233*a64d13efSStefano Zampini ierr = ISRestoreIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 234*a64d13efSStefano Zampini for (i=0, j=0; i<n_R; i++) { 235*a64d13efSStefano Zampini if (array_bool[idx_R_local[i]]) { 236*a64d13efSStefano Zampini aux_array1[j] = i; 237*a64d13efSStefano Zampini j++; 238*a64d13efSStefano Zampini } 239*a64d13efSStefano Zampini } 240*a64d13efSStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,j,aux_array1,PETSC_OWN_POINTER,&is_aux1);CHKERRQ(ierr); 241*a64d13efSStefano Zampini ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 242*a64d13efSStefano Zampini for (i=0, j=0; i<n_B; i++) { 243*a64d13efSStefano Zampini if (array_bool[is_indices[i]]) { 244*a64d13efSStefano Zampini aux_array2[j] = i; j++; 245*a64d13efSStefano Zampini } 246*a64d13efSStefano Zampini } 247*a64d13efSStefano Zampini ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 248*a64d13efSStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,j,aux_array2,PETSC_OWN_POINTER,&is_aux2);CHKERRQ(ierr); 249*a64d13efSStefano Zampini ierr = VecScatterCreate(pcbddc->vec1_R,is_aux1,pcis->vec1_B,is_aux2,&pcbddc->R_to_B);CHKERRQ(ierr); 250*a64d13efSStefano Zampini ierr = ISDestroy(&is_aux1);CHKERRQ(ierr); 251*a64d13efSStefano Zampini ierr = ISDestroy(&is_aux2);CHKERRQ(ierr); 252*a64d13efSStefano Zampini 253*a64d13efSStefano Zampini if (pcbddc->inexact_prec_type || pcbddc->dbg_flag ) { 254*a64d13efSStefano Zampini ierr = PetscMalloc(n_D*sizeof(PetscInt),&aux_array1);CHKERRQ(ierr); 255*a64d13efSStefano Zampini for (i=0, j=0; i<n_R; i++) { 256*a64d13efSStefano Zampini if (!array_bool[idx_R_local[i]]) { 257*a64d13efSStefano Zampini aux_array1[j] = i; 258*a64d13efSStefano Zampini j++; 259*a64d13efSStefano Zampini } 260*a64d13efSStefano Zampini } 261*a64d13efSStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,j,aux_array1,PETSC_OWN_POINTER,&is_aux1);CHKERRQ(ierr); 262*a64d13efSStefano Zampini ierr = VecScatterCreate(pcbddc->vec1_R,is_aux1,pcis->vec1_D,(IS)0,&pcbddc->R_to_D);CHKERRQ(ierr); 263*a64d13efSStefano Zampini ierr = ISDestroy(&is_aux1);CHKERRQ(ierr); 264*a64d13efSStefano Zampini } 265*a64d13efSStefano Zampini ierr = PetscFree(array_bool);CHKERRQ(ierr); 266*a64d13efSStefano Zampini *is_R_local_n = is_R_local; 267*a64d13efSStefano Zampini PetscFunctionReturn(0); 268*a64d13efSStefano Zampini } 269*a64d13efSStefano Zampini 270*a64d13efSStefano Zampini #undef __FUNCT__ 271304d26faSStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet" 272304d26faSStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool use) 273304d26faSStefano Zampini { 274304d26faSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 275304d26faSStefano Zampini 276304d26faSStefano Zampini PetscFunctionBegin; 277304d26faSStefano Zampini pcbddc->use_exact_dirichlet=use; 278304d26faSStefano Zampini PetscFunctionReturn(0); 279304d26faSStefano Zampini } 280304d26faSStefano Zampini 281304d26faSStefano Zampini #undef __FUNCT__ 282304d26faSStefano Zampini #define __FUNCT__ "PCBDDCSetUpLocalSolvers" 283304d26faSStefano Zampini PetscErrorCode PCBDDCSetUpLocalSolvers(PC pc, IS is_I_local, IS is_R_local) 284304d26faSStefano Zampini { 285304d26faSStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 286304d26faSStefano Zampini PC_IS *pcis = (PC_IS*)pc->data; 287304d26faSStefano Zampini PC pc_temp; 288304d26faSStefano Zampini Mat A_RR; 289304d26faSStefano Zampini Vec vec1,vec2,vec3; 290304d26faSStefano Zampini MatStructure matstruct; 291304d26faSStefano Zampini PetscScalar m_one = -1.0; 292304d26faSStefano Zampini PetscReal value; 293304d26faSStefano Zampini PetscInt n_D,n_R,use_exact,use_exact_reduced; 294304d26faSStefano Zampini PetscErrorCode ierr; 295304d26faSStefano Zampini 296304d26faSStefano Zampini PetscFunctionBegin; 297304d26faSStefano Zampini /* Creating PC contexts for local Dirichlet and Neumann problems */ 298304d26faSStefano Zampini ierr = PCGetOperators(pc,NULL,NULL,&matstruct);CHKERRQ(ierr); 299304d26faSStefano Zampini 300304d26faSStefano Zampini /* DIRICHLET PROBLEM */ 301ac78edfcSStefano Zampini /* Matrix for Dirichlet problem is pcis->A_II */ 302304d26faSStefano Zampini ierr = ISGetSize(is_I_local,&n_D);CHKERRQ(ierr); 303304d26faSStefano Zampini if (!pcbddc->ksp_D) { /* create object if not yet build */ 304304d26faSStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&pcbddc->ksp_D);CHKERRQ(ierr); 305304d26faSStefano Zampini ierr = PetscObjectIncrementTabLevel((PetscObject)pcbddc->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr); 306304d26faSStefano Zampini /* default */ 307304d26faSStefano Zampini ierr = KSPSetType(pcbddc->ksp_D,KSPPREONLY);CHKERRQ(ierr); 308304d26faSStefano Zampini ierr = KSPSetOptionsPrefix(pcbddc->ksp_D,"dirichlet_");CHKERRQ(ierr); 309304d26faSStefano Zampini ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr); 310304d26faSStefano Zampini ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr); 311304d26faSStefano Zampini ierr = PCFactorSetReuseFill(pc_temp,PETSC_TRUE);CHKERRQ(ierr); 312304d26faSStefano Zampini } 313304d26faSStefano Zampini ierr = KSPSetOperators(pcbddc->ksp_D,pcis->A_II,pcis->A_II,matstruct);CHKERRQ(ierr); 314304d26faSStefano Zampini /* Allow user's customization */ 315304d26faSStefano Zampini ierr = KSPSetFromOptions(pcbddc->ksp_D);CHKERRQ(ierr); 316304d26faSStefano Zampini /* umfpack interface has a bug when matrix dimension is zero. TODO solve from umfpack interface */ 317304d26faSStefano Zampini if (!n_D) { 318304d26faSStefano Zampini ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr); 319304d26faSStefano Zampini ierr = PCSetType(pc_temp,PCNONE);CHKERRQ(ierr); 320304d26faSStefano Zampini } 321304d26faSStefano Zampini /* Set Up KSP for Dirichlet problem of BDDC */ 322304d26faSStefano Zampini ierr = KSPSetUp(pcbddc->ksp_D);CHKERRQ(ierr); 323304d26faSStefano Zampini /* set ksp_D into pcis data */ 324304d26faSStefano Zampini ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr); 325304d26faSStefano Zampini ierr = PetscObjectReference((PetscObject)pcbddc->ksp_D);CHKERRQ(ierr); 326304d26faSStefano Zampini pcis->ksp_D = pcbddc->ksp_D; 327304d26faSStefano Zampini 328304d26faSStefano Zampini /* NEUMANN PROBLEM */ 329304d26faSStefano Zampini /* Matrix for Neumann problem is A_RR -> we need to create it */ 330304d26faSStefano Zampini ierr = ISGetSize(is_R_local,&n_R);CHKERRQ(ierr); 331304d26faSStefano Zampini ierr = MatGetSubMatrix(pcbddc->local_mat,is_R_local,is_R_local,MAT_INITIAL_MATRIX,&A_RR);CHKERRQ(ierr); 332304d26faSStefano Zampini if (!pcbddc->ksp_R) { /* create object if not yet build */ 333304d26faSStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&pcbddc->ksp_R);CHKERRQ(ierr); 334304d26faSStefano Zampini ierr = PetscObjectIncrementTabLevel((PetscObject)pcbddc->ksp_R,(PetscObject)pc,1);CHKERRQ(ierr); 335304d26faSStefano Zampini /* default */ 336304d26faSStefano Zampini ierr = KSPSetType(pcbddc->ksp_R,KSPPREONLY);CHKERRQ(ierr); 337304d26faSStefano Zampini ierr = KSPSetOptionsPrefix(pcbddc->ksp_R,"neumann_");CHKERRQ(ierr); 338304d26faSStefano Zampini ierr = KSPGetPC(pcbddc->ksp_R,&pc_temp);CHKERRQ(ierr); 339304d26faSStefano Zampini ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr); 340304d26faSStefano Zampini ierr = PCFactorSetReuseFill(pc_temp,PETSC_TRUE);CHKERRQ(ierr); 341304d26faSStefano Zampini } 342304d26faSStefano Zampini ierr = KSPSetOperators(pcbddc->ksp_R,A_RR,A_RR,matstruct);CHKERRQ(ierr); 343304d26faSStefano Zampini /* Allow user's customization */ 344304d26faSStefano Zampini ierr = KSPSetFromOptions(pcbddc->ksp_R);CHKERRQ(ierr); 345304d26faSStefano Zampini /* umfpack interface has a bug when matrix dimension is zero. TODO solve from umfpack interface */ 346304d26faSStefano Zampini if (!n_R) { 347304d26faSStefano Zampini ierr = KSPGetPC(pcbddc->ksp_R,&pc_temp);CHKERRQ(ierr); 348304d26faSStefano Zampini ierr = PCSetType(pc_temp,PCNONE);CHKERRQ(ierr); 349304d26faSStefano Zampini } 350304d26faSStefano Zampini /* Set Up KSP for Neumann problem of BDDC */ 351304d26faSStefano Zampini ierr = KSPSetUp(pcbddc->ksp_R);CHKERRQ(ierr); 352304d26faSStefano Zampini 353304d26faSStefano Zampini /* check Dirichlet and Neumann solvers and adapt them if a nullspace correction is needed */ 354304d26faSStefano Zampini 355304d26faSStefano Zampini /* Dirichlet */ 356304d26faSStefano Zampini ierr = MatGetVecs(pcis->A_II,&vec1,&vec2);CHKERRQ(ierr); 357304d26faSStefano Zampini ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr); 358304d26faSStefano Zampini ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr); 359304d26faSStefano Zampini ierr = MatMult(pcis->A_II,vec1,vec2);CHKERRQ(ierr); 360304d26faSStefano Zampini ierr = KSPSolve(pcbddc->ksp_D,vec2,vec3);CHKERRQ(ierr); 361304d26faSStefano Zampini ierr = VecAXPY(vec3,m_one,vec1);CHKERRQ(ierr); 362304d26faSStefano Zampini ierr = VecNorm(vec3,NORM_INFINITY,&value);CHKERRQ(ierr); 363304d26faSStefano Zampini ierr = VecDestroy(&vec1);CHKERRQ(ierr); 364304d26faSStefano Zampini ierr = VecDestroy(&vec2);CHKERRQ(ierr); 365304d26faSStefano Zampini ierr = VecDestroy(&vec3);CHKERRQ(ierr); 366304d26faSStefano Zampini /* need to be adapted? */ 367304d26faSStefano Zampini use_exact = (PetscAbsReal(value) > 1.e-4 ? 0 : 1); 368304d26faSStefano Zampini ierr = MPI_Allreduce(&use_exact,&use_exact_reduced,1,MPIU_INT,MPI_LAND,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 369304d26faSStefano Zampini ierr = PCBDDCSetUseExactDirichlet(pc,(PetscBool)use_exact_reduced);CHKERRQ(ierr); 370304d26faSStefano Zampini /* print info */ 371304d26faSStefano Zampini if (pcbddc->dbg_flag) { 372304d26faSStefano Zampini ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 373304d26faSStefano Zampini ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"--------------------------------------------------\n");CHKERRQ(ierr); 374304d26faSStefano Zampini ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Checking solution of Dirichlet and Neumann problems\n");CHKERRQ(ierr); 375304d26faSStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet solve = % 1.14e \n",PetscGlobalRank,value);CHKERRQ(ierr); 376304d26faSStefano Zampini ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 377304d26faSStefano Zampini } 378304d26faSStefano Zampini if (n_D && pcbddc->NullSpace && !use_exact_reduced && !pcbddc->inexact_prec_type) { 379304d26faSStefano Zampini ierr = PCBDDCNullSpaceAssembleCorrection(pc,is_I_local);CHKERRQ(ierr); 380304d26faSStefano Zampini } 381304d26faSStefano Zampini 382304d26faSStefano Zampini /* Neumann */ 383304d26faSStefano Zampini ierr = MatGetVecs(A_RR,&vec1,&vec2);CHKERRQ(ierr); 384304d26faSStefano Zampini ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr); 385304d26faSStefano Zampini ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr); 386304d26faSStefano Zampini ierr = MatMult(A_RR,vec1,vec2);CHKERRQ(ierr); 387304d26faSStefano Zampini ierr = KSPSolve(pcbddc->ksp_R,vec2,vec3);CHKERRQ(ierr); 388304d26faSStefano Zampini ierr = VecAXPY(vec3,m_one,vec1);CHKERRQ(ierr); 389304d26faSStefano Zampini ierr = VecNorm(vec3,NORM_INFINITY,&value);CHKERRQ(ierr); 390304d26faSStefano Zampini ierr = VecDestroy(&vec1);CHKERRQ(ierr); 391304d26faSStefano Zampini ierr = VecDestroy(&vec2);CHKERRQ(ierr); 392304d26faSStefano Zampini ierr = VecDestroy(&vec3);CHKERRQ(ierr); 393304d26faSStefano Zampini /* need to be adapted? */ 394304d26faSStefano Zampini use_exact = (PetscAbsReal(value) > 1.e-4 ? 0 : 1); 395304d26faSStefano Zampini if (PetscAbsReal(value) > 1.e-4) use_exact = 0; 396304d26faSStefano Zampini ierr = MPI_Allreduce(&use_exact,&use_exact_reduced,1,MPIU_INT,MPI_LAND,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 397304d26faSStefano Zampini /* print info */ 398304d26faSStefano Zampini if (pcbddc->dbg_flag) { 399304d26faSStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann solve = % 1.14e \n",PetscGlobalRank,value);CHKERRQ(ierr); 400304d26faSStefano Zampini ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 401304d26faSStefano Zampini } 402304d26faSStefano Zampini if (n_R && pcbddc->NullSpace && !use_exact_reduced) { /* is it the right logic? */ 403304d26faSStefano Zampini ierr = PCBDDCNullSpaceAssembleCorrection(pc,is_R_local);CHKERRQ(ierr); 404304d26faSStefano Zampini } 405304d26faSStefano Zampini 406304d26faSStefano Zampini /* free Neumann problem's matrix */ 407304d26faSStefano Zampini ierr = MatDestroy(&A_RR);CHKERRQ(ierr); 408304d26faSStefano Zampini PetscFunctionReturn(0); 409304d26faSStefano Zampini } 410304d26faSStefano Zampini 411304d26faSStefano Zampini #undef __FUNCT__ 412674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSolveSaddlePoint" 413674ae819SStefano Zampini static PetscErrorCode PCBDDCSolveSaddlePoint(PC pc) 414674ae819SStefano Zampini { 415674ae819SStefano Zampini PetscErrorCode ierr; 416674ae819SStefano Zampini PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 417674ae819SStefano Zampini 418674ae819SStefano Zampini PetscFunctionBegin; 419674ae819SStefano Zampini ierr = KSPSolve(pcbddc->ksp_R,pcbddc->vec1_R,pcbddc->vec2_R);CHKERRQ(ierr); 420674ae819SStefano Zampini if (pcbddc->local_auxmat1) { 421674ae819SStefano Zampini ierr = MatMult(pcbddc->local_auxmat1,pcbddc->vec2_R,pcbddc->vec1_C);CHKERRQ(ierr); 422674ae819SStefano Zampini ierr = MatMultAdd(pcbddc->local_auxmat2,pcbddc->vec1_C,pcbddc->vec2_R,pcbddc->vec2_R);CHKERRQ(ierr); 423674ae819SStefano Zampini } 424674ae819SStefano Zampini PetscFunctionReturn(0); 425674ae819SStefano Zampini } 426674ae819SStefano Zampini 427674ae819SStefano Zampini #undef __FUNCT__ 428674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplyInterfacePreconditioner" 429674ae819SStefano Zampini PetscErrorCode PCBDDCApplyInterfacePreconditioner(PC pc) 430674ae819SStefano Zampini { 431674ae819SStefano Zampini PetscErrorCode ierr; 432674ae819SStefano Zampini PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 433674ae819SStefano Zampini PC_IS* pcis = (PC_IS*) (pc->data); 434674ae819SStefano Zampini const PetscScalar zero = 0.0; 435674ae819SStefano Zampini 436674ae819SStefano Zampini PetscFunctionBegin; 43715aaf578SStefano Zampini /* Application of PHI^T (or PSI^T) */ 43815aaf578SStefano Zampini if (pcbddc->coarse_psi_B) { 43915aaf578SStefano Zampini ierr = MatMultTranspose(pcbddc->coarse_psi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr); 44015aaf578SStefano Zampini if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_psi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); } 44115aaf578SStefano Zampini } else { 442674ae819SStefano Zampini ierr = MatMultTranspose(pcbddc->coarse_phi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr); 443674ae819SStefano Zampini if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_phi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); } 44415aaf578SStefano Zampini } 445674ae819SStefano Zampini /* Scatter data of coarse_rhs */ 446674ae819SStefano Zampini if (pcbddc->coarse_rhs) { ierr = VecSet(pcbddc->coarse_rhs,zero);CHKERRQ(ierr); } 447674ae819SStefano Zampini ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 448674ae819SStefano Zampini 449674ae819SStefano Zampini /* Local solution on R nodes */ 450674ae819SStefano Zampini ierr = VecSet(pcbddc->vec1_R,zero);CHKERRQ(ierr); 451674ae819SStefano Zampini ierr = VecScatterBegin(pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 452674ae819SStefano Zampini ierr = VecScatterEnd (pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 453674ae819SStefano Zampini if (pcbddc->inexact_prec_type) { 454674ae819SStefano Zampini ierr = VecScatterBegin(pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 455674ae819SStefano Zampini ierr = VecScatterEnd (pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 456674ae819SStefano Zampini } 457674ae819SStefano Zampini ierr = PCBDDCSolveSaddlePoint(pc);CHKERRQ(ierr); 458674ae819SStefano Zampini ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr); 459674ae819SStefano Zampini ierr = VecScatterBegin(pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 460674ae819SStefano Zampini ierr = VecScatterEnd (pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 461674ae819SStefano Zampini if (pcbddc->inexact_prec_type) { 462674ae819SStefano Zampini ierr = VecScatterBegin(pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 463674ae819SStefano Zampini ierr = VecScatterEnd (pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 464674ae819SStefano Zampini } 465674ae819SStefano Zampini 466674ae819SStefano Zampini /* Coarse solution */ 467674ae819SStefano Zampini ierr = PCBDDCScatterCoarseDataEnd(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 468674ae819SStefano Zampini if (pcbddc->coarse_rhs) { /* TODO remove null space when doing multilevel */ 469674ae819SStefano Zampini ierr = KSPSolve(pcbddc->coarse_ksp,pcbddc->coarse_rhs,pcbddc->coarse_vec);CHKERRQ(ierr); 470674ae819SStefano Zampini } 471674ae819SStefano Zampini ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 472674ae819SStefano Zampini ierr = PCBDDCScatterCoarseDataEnd (pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 473674ae819SStefano Zampini 474674ae819SStefano Zampini /* Sum contributions from two levels */ 475674ae819SStefano Zampini ierr = MatMultAdd(pcbddc->coarse_phi_B,pcbddc->vec1_P,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr); 476674ae819SStefano Zampini if (pcbddc->inexact_prec_type) { ierr = MatMultAdd(pcbddc->coarse_phi_D,pcbddc->vec1_P,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 477674ae819SStefano Zampini PetscFunctionReturn(0); 478674ae819SStefano Zampini } 479674ae819SStefano Zampini 480674ae819SStefano Zampini #undef __FUNCT__ 481674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScatterCoarseDataBegin" 482674ae819SStefano Zampini PetscErrorCode PCBDDCScatterCoarseDataBegin(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode) 483674ae819SStefano Zampini { 484674ae819SStefano Zampini PetscErrorCode ierr; 485674ae819SStefano Zampini PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 486674ae819SStefano Zampini 487674ae819SStefano Zampini PetscFunctionBegin; 488674ae819SStefano Zampini switch (pcbddc->coarse_communications_type) { 489674ae819SStefano Zampini case SCATTERS_BDDC: 490674ae819SStefano Zampini ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr); 491674ae819SStefano Zampini break; 492674ae819SStefano Zampini case GATHERS_BDDC: 493674ae819SStefano Zampini break; 494674ae819SStefano Zampini } 495674ae819SStefano Zampini PetscFunctionReturn(0); 496674ae819SStefano Zampini } 497674ae819SStefano Zampini 498674ae819SStefano Zampini #undef __FUNCT__ 499674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScatterCoarseDataEnd" 500674ae819SStefano Zampini PetscErrorCode PCBDDCScatterCoarseDataEnd(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode) 501674ae819SStefano Zampini { 502674ae819SStefano Zampini PetscErrorCode ierr; 503674ae819SStefano Zampini PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 504674ae819SStefano Zampini PetscScalar* array_to; 505674ae819SStefano Zampini PetscScalar* array_from; 506674ae819SStefano Zampini MPI_Comm comm; 507674ae819SStefano Zampini PetscInt i; 508674ae819SStefano Zampini 509674ae819SStefano Zampini PetscFunctionBegin; 510674ae819SStefano Zampini ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 511674ae819SStefano Zampini switch (pcbddc->coarse_communications_type) { 512674ae819SStefano Zampini case SCATTERS_BDDC: 513674ae819SStefano Zampini ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr); 514674ae819SStefano Zampini break; 515674ae819SStefano Zampini case GATHERS_BDDC: 516674ae819SStefano Zampini if (vec_from) { 517674ae819SStefano Zampini ierr = VecGetArray(vec_from,&array_from);CHKERRQ(ierr); 518674ae819SStefano Zampini } 519674ae819SStefano Zampini if (vec_to) { 520674ae819SStefano Zampini ierr = VecGetArray(vec_to,&array_to);CHKERRQ(ierr); 521674ae819SStefano Zampini } 522674ae819SStefano Zampini switch(pcbddc->coarse_problem_type){ 523674ae819SStefano Zampini case SEQUENTIAL_BDDC: 524674ae819SStefano Zampini if (smode == SCATTER_FORWARD) { 525674ae819SStefano Zampini ierr = MPI_Gatherv(&array_from[0],pcbddc->local_primal_size,MPIU_SCALAR,&pcbddc->replicated_local_primal_values[0],pcbddc->local_primal_sizes,pcbddc->local_primal_displacements,MPIU_SCALAR,0,comm);CHKERRQ(ierr); 526674ae819SStefano Zampini if (vec_to) { 527674ae819SStefano Zampini if (imode == ADD_VALUES) { 528674ae819SStefano Zampini for (i=0;i<pcbddc->replicated_primal_size;i++) { 529674ae819SStefano Zampini array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i]; 530674ae819SStefano Zampini } 531674ae819SStefano Zampini } else { 532674ae819SStefano Zampini for (i=0;i<pcbddc->replicated_primal_size;i++) { 533674ae819SStefano Zampini array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i]; 534674ae819SStefano Zampini } 535674ae819SStefano Zampini } 536674ae819SStefano Zampini } 537674ae819SStefano Zampini } else { 538674ae819SStefano Zampini if (vec_from) { 539674ae819SStefano Zampini if (imode == ADD_VALUES) { 540674ae819SStefano Zampini MPI_Comm vec_from_comm; 541674ae819SStefano Zampini ierr = PetscObjectGetComm((PetscObject)(vec_from),&vec_from_comm);CHKERRQ(ierr); 542674ae819SStefano Zampini SETERRQ2(vec_from_comm,PETSC_ERR_SUP,"Unsupported insert mode ADD_VALUES for SCATTER_REVERSE in %s for case %d\n",__FUNCT__,pcbddc->coarse_problem_type); 543674ae819SStefano Zampini } 544674ae819SStefano Zampini for (i=0;i<pcbddc->replicated_primal_size;i++) { 545674ae819SStefano Zampini pcbddc->replicated_local_primal_values[i]=array_from[pcbddc->replicated_local_primal_indices[i]]; 546674ae819SStefano Zampini } 547674ae819SStefano Zampini } 548674ae819SStefano Zampini ierr = MPI_Scatterv(&pcbddc->replicated_local_primal_values[0],pcbddc->local_primal_sizes,pcbddc->local_primal_displacements,MPIU_SCALAR,&array_to[0],pcbddc->local_primal_size,MPIU_SCALAR,0,comm);CHKERRQ(ierr); 549674ae819SStefano Zampini } 550674ae819SStefano Zampini break; 551674ae819SStefano Zampini case REPLICATED_BDDC: 552674ae819SStefano Zampini if (smode == SCATTER_FORWARD) { 553674ae819SStefano Zampini ierr = MPI_Allgatherv(&array_from[0],pcbddc->local_primal_size,MPIU_SCALAR,&pcbddc->replicated_local_primal_values[0],pcbddc->local_primal_sizes,pcbddc->local_primal_displacements,MPIU_SCALAR,comm);CHKERRQ(ierr); 554674ae819SStefano Zampini if (imode == ADD_VALUES) { 555674ae819SStefano Zampini for (i=0;i<pcbddc->replicated_primal_size;i++) { 556674ae819SStefano Zampini array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i]; 557674ae819SStefano Zampini } 558674ae819SStefano Zampini } else { 559674ae819SStefano Zampini for (i=0;i<pcbddc->replicated_primal_size;i++) { 560674ae819SStefano Zampini array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i]; 561674ae819SStefano Zampini } 562674ae819SStefano Zampini } 563674ae819SStefano Zampini } else { /* no communications needed for SCATTER_REVERSE since needed data is already present */ 564674ae819SStefano Zampini if (imode == ADD_VALUES) { 565674ae819SStefano Zampini for (i=0;i<pcbddc->local_primal_size;i++) { 566674ae819SStefano Zampini array_to[i]+=array_from[pcbddc->local_primal_indices[i]]; 567674ae819SStefano Zampini } 568674ae819SStefano Zampini } else { 569674ae819SStefano Zampini for (i=0;i<pcbddc->local_primal_size;i++) { 570674ae819SStefano Zampini array_to[i]=array_from[pcbddc->local_primal_indices[i]]; 571674ae819SStefano Zampini } 572674ae819SStefano Zampini } 573674ae819SStefano Zampini } 574674ae819SStefano Zampini break; 575674ae819SStefano Zampini case MULTILEVEL_BDDC: 576674ae819SStefano Zampini break; 577674ae819SStefano Zampini case PARALLEL_BDDC: 578674ae819SStefano Zampini break; 579674ae819SStefano Zampini } 580674ae819SStefano Zampini if (vec_from) { 581674ae819SStefano Zampini ierr = VecRestoreArray(vec_from,&array_from);CHKERRQ(ierr); 582674ae819SStefano Zampini } 583674ae819SStefano Zampini if (vec_to) { 584674ae819SStefano Zampini ierr = VecRestoreArray(vec_to,&array_to);CHKERRQ(ierr); 585674ae819SStefano Zampini } 586674ae819SStefano Zampini break; 587674ae819SStefano Zampini } 588674ae819SStefano Zampini PetscFunctionReturn(0); 589674ae819SStefano Zampini } 590674ae819SStefano Zampini 591674ae819SStefano Zampini #undef __FUNCT__ 592674ae819SStefano Zampini #define __FUNCT__ "PCBDDCConstraintsSetUp" 593674ae819SStefano Zampini PetscErrorCode PCBDDCConstraintsSetUp(PC pc) 594674ae819SStefano Zampini { 595674ae819SStefano Zampini PetscErrorCode ierr; 596674ae819SStefano Zampini PC_IS* pcis = (PC_IS*)(pc->data); 597674ae819SStefano Zampini PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 598674ae819SStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 599674ae819SStefano Zampini PetscInt *nnz,*is_indices; 600674ae819SStefano Zampini PetscScalar *temp_quadrature_constraint; 601674ae819SStefano Zampini PetscInt *temp_indices,*temp_indices_to_constraint,*temp_indices_to_constraint_B,*local_to_B; 602674ae819SStefano Zampini PetscInt local_primal_size,i,j,k,total_counts,max_size_of_constraint; 603674ae819SStefano Zampini PetscInt n_vertices,size_of_constraint; 6045b08dc53SStefano Zampini PetscReal real_value; 605674ae819SStefano Zampini PetscBool nnsp_has_cnst=PETSC_FALSE,use_nnsp_true=pcbddc->use_nnsp_true; 606674ae819SStefano Zampini PetscInt nnsp_size=0,nnsp_addone=0,temp_constraints,temp_start_ptr,n_ISForFaces,n_ISForEdges; 607674ae819SStefano Zampini IS *used_IS,ISForVertices,*ISForFaces,*ISForEdges; 608674ae819SStefano Zampini MatType impMatType=MATSEQAIJ; 609674ae819SStefano Zampini PetscBLASInt Bs,Bt,lwork,lierr; 610674ae819SStefano Zampini PetscReal tol=1.0e-8; 611674ae819SStefano Zampini MatNullSpace nearnullsp; 612674ae819SStefano Zampini const Vec *nearnullvecs; 613674ae819SStefano Zampini Vec *localnearnullsp; 614674ae819SStefano Zampini PetscScalar *work,*temp_basis,*array_vector,*correlation_mat; 615674ae819SStefano Zampini PetscReal *rwork,*singular_vals; 616674ae819SStefano Zampini PetscBLASInt Bone=1,*ipiv; 617674ae819SStefano Zampini Vec temp_vec; 618674ae819SStefano Zampini Mat temp_mat; 619674ae819SStefano Zampini KSP temp_ksp; 620674ae819SStefano Zampini PC temp_pc; 621674ae819SStefano Zampini PetscInt s,start_constraint,dual_dofs; 622674ae819SStefano Zampini PetscBool compute_submatrix,useksp=PETSC_FALSE; 623674ae819SStefano Zampini PetscInt *aux_primal_permutation,*aux_primal_numbering; 62451b0f6cfSStefano Zampini PetscBool boolforchange,*change_basis; 625674ae819SStefano Zampini /* some ugly conditional declarations */ 626674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD) 627674ae819SStefano Zampini PetscScalar one=1.0,zero=0.0; 628674ae819SStefano Zampini PetscInt ii; 629674ae819SStefano Zampini PetscScalar *singular_vectors; 630674ae819SStefano Zampini PetscBLASInt *iwork,*ifail; 631674ae819SStefano Zampini PetscReal dummy_real,abs_tol; 632674ae819SStefano Zampini PetscBLASInt eigs_found; 633674ae819SStefano Zampini #endif 634674ae819SStefano Zampini PetscBLASInt dummy_int; 635674ae819SStefano Zampini PetscScalar dummy_scalar; 636674ae819SStefano Zampini PetscBool used_vertex,get_faces,get_edges,get_vertices; 637674ae819SStefano Zampini 638674ae819SStefano Zampini PetscFunctionBegin; 639674ae819SStefano Zampini /* Get index sets for faces, edges and vertices from graph */ 640674ae819SStefano Zampini get_faces = PETSC_TRUE; 641674ae819SStefano Zampini get_edges = PETSC_TRUE; 642674ae819SStefano Zampini get_vertices = PETSC_TRUE; 643674ae819SStefano Zampini if (pcbddc->vertices_flag) { 644674ae819SStefano Zampini get_faces = PETSC_FALSE; 645674ae819SStefano Zampini get_edges = PETSC_FALSE; 646674ae819SStefano Zampini } 647674ae819SStefano Zampini if (pcbddc->constraints_flag) { 648674ae819SStefano Zampini get_vertices = PETSC_FALSE; 649674ae819SStefano Zampini } 650674ae819SStefano Zampini if (pcbddc->faces_flag) { 651674ae819SStefano Zampini get_edges = PETSC_FALSE; 652674ae819SStefano Zampini } 653674ae819SStefano Zampini if (pcbddc->edges_flag) { 654674ae819SStefano Zampini get_faces = PETSC_FALSE; 655674ae819SStefano Zampini } 656674ae819SStefano Zampini /* default */ 657674ae819SStefano Zampini if (!get_faces && !get_edges && !get_vertices) { 658674ae819SStefano Zampini get_vertices = PETSC_TRUE; 659674ae819SStefano Zampini } 660674ae819SStefano Zampini ierr = PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,get_faces,get_edges,get_vertices,&n_ISForFaces,&ISForFaces,&n_ISForEdges,&ISForEdges,&ISForVertices); 661674ae819SStefano Zampini if (pcbddc->dbg_flag) { 662674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"--------------------------------------------------------------\n");CHKERRQ(ierr); 663674ae819SStefano Zampini i = 0; 664674ae819SStefano Zampini if (ISForVertices) { 665674ae819SStefano Zampini ierr = ISGetSize(ISForVertices,&i);CHKERRQ(ierr); 666674ae819SStefano Zampini } 667674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate vertices\n",PetscGlobalRank,i);CHKERRQ(ierr); 668674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate edges\n",PetscGlobalRank,n_ISForEdges);CHKERRQ(ierr); 66915aaf578SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate faces\n",PetscGlobalRank,n_ISForFaces);CHKERRQ(ierr); 670674ae819SStefano Zampini ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 671674ae819SStefano Zampini } 672674ae819SStefano Zampini /* check if near null space is attached to global mat */ 673674ae819SStefano Zampini ierr = MatGetNearNullSpace(pc->pmat,&nearnullsp);CHKERRQ(ierr); 674674ae819SStefano Zampini if (nearnullsp) { 675674ae819SStefano Zampini ierr = MatNullSpaceGetVecs(nearnullsp,&nnsp_has_cnst,&nnsp_size,&nearnullvecs);CHKERRQ(ierr); 676674ae819SStefano Zampini } else { /* if near null space is not provided it uses constants */ 677674ae819SStefano Zampini nnsp_has_cnst = PETSC_TRUE; 678674ae819SStefano Zampini use_nnsp_true = PETSC_TRUE; 679674ae819SStefano Zampini } 680674ae819SStefano Zampini if (nnsp_has_cnst) { 681674ae819SStefano Zampini nnsp_addone = 1; 682674ae819SStefano Zampini } 683674ae819SStefano Zampini /* 684674ae819SStefano Zampini Evaluate maximum storage size needed by the procedure 685674ae819SStefano Zampini - temp_indices will contain start index of each constraint stored as follows 686674ae819SStefano Zampini - temp_indices_to_constraint [temp_indices[i],...,temp[indices[i+1]-1] will contain the indices (in local numbering) on which the constraint acts 687674ae819SStefano Zampini - temp_indices_to_constraint_B[temp_indices[i],...,temp[indices[i+1]-1] will contain the indices (in boundary numbering) on which the constraint acts 688674ae819SStefano Zampini - temp_quadrature_constraint [temp_indices[i],...,temp[indices[i+1]-1] will contain the scalars representing the constraint itself 689674ae819SStefano Zampini */ 690674ae819SStefano Zampini total_counts = n_ISForFaces+n_ISForEdges; 691674ae819SStefano Zampini total_counts *= (nnsp_addone+nnsp_size); 692674ae819SStefano Zampini n_vertices = 0; 693674ae819SStefano Zampini if (ISForVertices) { 694674ae819SStefano Zampini ierr = ISGetSize(ISForVertices,&n_vertices);CHKERRQ(ierr); 695674ae819SStefano Zampini } 696674ae819SStefano Zampini total_counts += n_vertices; 697674ae819SStefano Zampini ierr = PetscMalloc((total_counts+1)*sizeof(PetscInt),&temp_indices);CHKERRQ(ierr); 698674ae819SStefano Zampini ierr = PetscMalloc((total_counts+1)*sizeof(PetscBool),&change_basis);CHKERRQ(ierr); 699674ae819SStefano Zampini total_counts = 0; 700674ae819SStefano Zampini max_size_of_constraint = 0; 701674ae819SStefano Zampini for (i=0;i<n_ISForEdges+n_ISForFaces;i++) { 702674ae819SStefano Zampini if (i<n_ISForEdges) { 703674ae819SStefano Zampini used_IS = &ISForEdges[i]; 704674ae819SStefano Zampini } else { 705674ae819SStefano Zampini used_IS = &ISForFaces[i-n_ISForEdges]; 706674ae819SStefano Zampini } 707674ae819SStefano Zampini ierr = ISGetSize(*used_IS,&j);CHKERRQ(ierr); 708674ae819SStefano Zampini total_counts += j; 709674ae819SStefano Zampini max_size_of_constraint = PetscMax(j,max_size_of_constraint); 710674ae819SStefano Zampini } 711674ae819SStefano Zampini total_counts *= (nnsp_addone+nnsp_size); 712674ae819SStefano Zampini total_counts += n_vertices; 713674ae819SStefano Zampini ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&temp_quadrature_constraint);CHKERRQ(ierr); 714674ae819SStefano Zampini ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint);CHKERRQ(ierr); 715674ae819SStefano Zampini ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint_B);CHKERRQ(ierr); 716674ae819SStefano Zampini ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&local_to_B);CHKERRQ(ierr); 717674ae819SStefano Zampini ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 718674ae819SStefano Zampini for (i=0;i<pcis->n;i++) { 719674ae819SStefano Zampini local_to_B[i]=-1; 720674ae819SStefano Zampini } 721674ae819SStefano Zampini for (i=0;i<pcis->n_B;i++) { 722674ae819SStefano Zampini local_to_B[is_indices[i]]=i; 723674ae819SStefano Zampini } 724674ae819SStefano Zampini ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 725674ae819SStefano Zampini 726674ae819SStefano Zampini /* First we issue queries to allocate optimal workspace for LAPACKgesvd or LAPACKsyev/LAPACKheev */ 727674ae819SStefano Zampini rwork = 0; 728674ae819SStefano Zampini work = 0; 729674ae819SStefano Zampini singular_vals = 0; 730674ae819SStefano Zampini temp_basis = 0; 731674ae819SStefano Zampini correlation_mat = 0; 732674ae819SStefano Zampini if (!pcbddc->use_nnsp_true) { 733674ae819SStefano Zampini PetscScalar temp_work; 734674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD) 735674ae819SStefano Zampini /* POD */ 736674ae819SStefano Zampini PetscInt max_n; 737674ae819SStefano Zampini max_n = nnsp_addone+nnsp_size; 738674ae819SStefano Zampini /* using some techniques borrowed from Proper Orthogonal Decomposition */ 739674ae819SStefano Zampini ierr = PetscMalloc(max_n*max_n*sizeof(PetscScalar),&correlation_mat);CHKERRQ(ierr); 740674ae819SStefano Zampini ierr = PetscMalloc(max_n*max_n*sizeof(PetscScalar),&singular_vectors);CHKERRQ(ierr); 741674ae819SStefano Zampini ierr = PetscMalloc(max_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr); 742674ae819SStefano Zampini ierr = PetscMalloc(max_size_of_constraint*(nnsp_addone+nnsp_size)*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr); 743674ae819SStefano Zampini #if defined(PETSC_USE_COMPLEX) 744674ae819SStefano Zampini ierr = PetscMalloc(3*max_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr); 745674ae819SStefano Zampini #endif 746674ae819SStefano Zampini ierr = PetscMalloc(5*max_n*sizeof(PetscBLASInt),&iwork);CHKERRQ(ierr); 747674ae819SStefano Zampini ierr = PetscMalloc(max_n*sizeof(PetscBLASInt),&ifail);CHKERRQ(ierr); 748674ae819SStefano Zampini /* now we evaluate the optimal workspace using query with lwork=-1 */ 749674ae819SStefano Zampini ierr = PetscBLASIntCast(max_n,&Bt);CHKERRQ(ierr); 750674ae819SStefano Zampini lwork=-1; 751674ae819SStefano Zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 752674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX) 753674ae819SStefano Zampini abs_tol=1.e-8; 754674ae819SStefano Zampini /* LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,&lierr); */ 755f7c40c41SStefano Zampini PetscStackCallBLAS("LAPACKsyevx",LAPACKsyevx_("V","A","U",&Bt,correlation_mat,&Bt,&dummy_real,&dummy_real,&dummy_int,&dummy_int,&abs_tol,&eigs_found,singular_vals,singular_vectors,&Bt,&temp_work,&lwork,iwork,ifail,&lierr)); 756674ae819SStefano Zampini #else 757674ae819SStefano Zampini /* LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,rwork,&lierr); */ 758674ae819SStefano Zampini /* LAPACK call is missing here! TODO */ 759674ae819SStefano Zampini SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for complexes when PETSC_MISSING_GESVD = 1"); 760674ae819SStefano Zampini #endif 761674ae819SStefano Zampini if ( lierr ) { 762674ae819SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYEVX Lapack routine %d",(int)lierr); 763674ae819SStefano Zampini } 764674ae819SStefano Zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 765674ae819SStefano Zampini #else /* on missing GESVD */ 766674ae819SStefano Zampini /* SVD */ 767674ae819SStefano Zampini PetscInt max_n,min_n; 768674ae819SStefano Zampini max_n = max_size_of_constraint; 769674ae819SStefano Zampini min_n = nnsp_addone+nnsp_size; 770674ae819SStefano Zampini if (max_size_of_constraint < ( nnsp_addone+nnsp_size ) ) { 771674ae819SStefano Zampini min_n = max_size_of_constraint; 772674ae819SStefano Zampini max_n = nnsp_addone+nnsp_size; 773674ae819SStefano Zampini } 774674ae819SStefano Zampini ierr = PetscMalloc(min_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr); 775674ae819SStefano Zampini #if defined(PETSC_USE_COMPLEX) 776674ae819SStefano Zampini ierr = PetscMalloc(5*min_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr); 777674ae819SStefano Zampini #endif 778674ae819SStefano Zampini /* now we evaluate the optimal workspace using query with lwork=-1 */ 779674ae819SStefano Zampini lwork=-1; 780674ae819SStefano Zampini ierr = PetscBLASIntCast(max_n,&Bs);CHKERRQ(ierr); 781674ae819SStefano Zampini ierr = PetscBLASIntCast(min_n,&Bt);CHKERRQ(ierr); 782674ae819SStefano Zampini dummy_int = Bs; 783674ae819SStefano Zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 784674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX) 785f7c40c41SStefano Zampini PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[0],&Bs,singular_vals,&dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,&temp_work,&lwork,&lierr)); 786674ae819SStefano Zampini #else 787f7c40c41SStefano Zampini PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[0],&Bs,singular_vals,&dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,&temp_work,&lwork,rwork,&lierr)); 788674ae819SStefano Zampini #endif 789674ae819SStefano Zampini if ( lierr ) { 790674ae819SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SVD Lapack routine %d",(int)lierr); 791674ae819SStefano Zampini } 792674ae819SStefano Zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 793674ae819SStefano Zampini #endif 794674ae819SStefano Zampini /* Allocate optimal workspace */ 795674ae819SStefano Zampini ierr = PetscBLASIntCast((PetscInt)PetscRealPart(temp_work),&lwork);CHKERRQ(ierr); 796674ae819SStefano Zampini total_counts = (PetscInt)lwork; 797674ae819SStefano Zampini ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&work);CHKERRQ(ierr); 798674ae819SStefano Zampini } 799674ae819SStefano Zampini /* get local part of global near null space vectors */ 800674ae819SStefano Zampini ierr = PetscMalloc(nnsp_size*sizeof(Vec),&localnearnullsp);CHKERRQ(ierr); 801674ae819SStefano Zampini for (k=0;k<nnsp_size;k++) { 802674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_N,&localnearnullsp[k]);CHKERRQ(ierr); 803674ae819SStefano Zampini ierr = VecScatterBegin(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 804674ae819SStefano Zampini ierr = VecScatterEnd(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 805674ae819SStefano Zampini } 806674ae819SStefano Zampini /* Now we can loop on constraining sets */ 807674ae819SStefano Zampini total_counts = 0; 808674ae819SStefano Zampini temp_indices[0] = 0; 809674ae819SStefano Zampini /* vertices */ 810674ae819SStefano Zampini if (ISForVertices) { 811674ae819SStefano Zampini ierr = ISGetIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr); 812674ae819SStefano Zampini if (nnsp_has_cnst) { /* consider all vertices */ 813674ae819SStefano Zampini for (i=0;i<n_vertices;i++) { 814674ae819SStefano Zampini temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i]; 815674ae819SStefano Zampini temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]]; 816674ae819SStefano Zampini temp_quadrature_constraint[temp_indices[total_counts]]=1.0; 817674ae819SStefano Zampini temp_indices[total_counts+1]=temp_indices[total_counts]+1; 818674ae819SStefano Zampini change_basis[total_counts]=PETSC_FALSE; 819674ae819SStefano Zampini total_counts++; 820674ae819SStefano Zampini } 821674ae819SStefano Zampini } else { /* consider vertices for which exist at least a localnearnullsp which is not null there */ 822674ae819SStefano Zampini for (i=0;i<n_vertices;i++) { 823674ae819SStefano Zampini used_vertex=PETSC_FALSE; 824674ae819SStefano Zampini k=0; 825674ae819SStefano Zampini while (!used_vertex && k<nnsp_size) { 826674ae819SStefano Zampini ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr); 827674ae819SStefano Zampini if (PetscAbsScalar(array_vector[is_indices[i]])>0.0) { 828674ae819SStefano Zampini temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i]; 829674ae819SStefano Zampini temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]]; 830674ae819SStefano Zampini temp_quadrature_constraint[temp_indices[total_counts]]=1.0; 831674ae819SStefano Zampini temp_indices[total_counts+1]=temp_indices[total_counts]+1; 832674ae819SStefano Zampini change_basis[total_counts]=PETSC_FALSE; 833674ae819SStefano Zampini total_counts++; 834674ae819SStefano Zampini used_vertex=PETSC_TRUE; 835674ae819SStefano Zampini } 836674ae819SStefano Zampini ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr); 837674ae819SStefano Zampini k++; 838674ae819SStefano Zampini } 839674ae819SStefano Zampini } 840674ae819SStefano Zampini } 841674ae819SStefano Zampini ierr = ISRestoreIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr); 842674ae819SStefano Zampini n_vertices = total_counts; 843674ae819SStefano Zampini } 844674ae819SStefano Zampini /* edges and faces */ 845674ae819SStefano Zampini for (i=0;i<n_ISForEdges+n_ISForFaces;i++) { 846674ae819SStefano Zampini if (i<n_ISForEdges) { 847674ae819SStefano Zampini used_IS = &ISForEdges[i]; 84851b0f6cfSStefano Zampini boolforchange = pcbddc->use_change_of_basis; 849674ae819SStefano Zampini } else { 850674ae819SStefano Zampini used_IS = &ISForFaces[i-n_ISForEdges]; 85151b0f6cfSStefano Zampini boolforchange = pcbddc->use_change_on_faces; 852674ae819SStefano Zampini } 853674ae819SStefano Zampini temp_constraints = 0; /* zero the number of constraints I have on this conn comp */ 854674ae819SStefano Zampini temp_start_ptr = total_counts; /* need to know the starting index of constraints stored */ 855674ae819SStefano Zampini ierr = ISGetSize(*used_IS,&size_of_constraint);CHKERRQ(ierr); 856674ae819SStefano Zampini ierr = ISGetIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 85751b0f6cfSStefano Zampini /* HACK: change of basis should not performed on local periodic nodes */ 85851b0f6cfSStefano Zampini if (pcbddc->mat_graph->mirrors && pcbddc->mat_graph->mirrors[is_indices[0]]) { 85951b0f6cfSStefano Zampini boolforchange = PETSC_FALSE; 86051b0f6cfSStefano Zampini } 861674ae819SStefano Zampini if (nnsp_has_cnst) { 8625b08dc53SStefano Zampini PetscScalar quad_value; 863674ae819SStefano Zampini temp_constraints++; 864674ae819SStefano Zampini quad_value = (PetscScalar)(1.0/PetscSqrtReal((PetscReal)size_of_constraint)); 865674ae819SStefano Zampini for (j=0;j<size_of_constraint;j++) { 866674ae819SStefano Zampini temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j]; 867674ae819SStefano Zampini temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]]; 868674ae819SStefano Zampini temp_quadrature_constraint[temp_indices[total_counts]+j]=quad_value; 869674ae819SStefano Zampini } 870674ae819SStefano Zampini temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint; /* store new starting point */ 87151b0f6cfSStefano Zampini change_basis[total_counts]=boolforchange; 872674ae819SStefano Zampini total_counts++; 873674ae819SStefano Zampini } 874674ae819SStefano Zampini for (k=0;k<nnsp_size;k++) { 875674ae819SStefano Zampini ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr); 876674ae819SStefano Zampini for (j=0;j<size_of_constraint;j++) { 877674ae819SStefano Zampini temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j]; 878674ae819SStefano Zampini temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]]; 879674ae819SStefano Zampini temp_quadrature_constraint[temp_indices[total_counts]+j]=array_vector[is_indices[j]]; 880674ae819SStefano Zampini } 881674ae819SStefano Zampini ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr); 8825b08dc53SStefano Zampini real_value = 1.0; 883674ae819SStefano Zampini if (use_nnsp_true) { /* check if array is null on the connected component in case use_nnsp_true has been requested */ 884674ae819SStefano Zampini ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr); 8855b08dc53SStefano Zampini PetscStackCallBLAS("BLASasum",real_value = BLASasum_(&Bs,&temp_quadrature_constraint[temp_indices[total_counts]],&Bone)); 886674ae819SStefano Zampini } 8875b08dc53SStefano Zampini if (real_value > 0.0) { /* keep indices and values */ 888674ae819SStefano Zampini temp_constraints++; 889674ae819SStefano Zampini temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint; /* store new starting point */ 89051b0f6cfSStefano Zampini change_basis[total_counts]=boolforchange; 891674ae819SStefano Zampini total_counts++; 892674ae819SStefano Zampini } 893674ae819SStefano Zampini } 894674ae819SStefano Zampini ierr = ISRestoreIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 895674ae819SStefano Zampini /* perform SVD on the constraint if use_nnsp_true has not be requested by the user */ 896674ae819SStefano Zampini if (!use_nnsp_true) { 897674ae819SStefano Zampini ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr); 898674ae819SStefano Zampini ierr = PetscBLASIntCast(temp_constraints,&Bt);CHKERRQ(ierr); 899674ae819SStefano Zampini 900674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD) 901674ae819SStefano Zampini ierr = PetscMemzero(correlation_mat,Bt*Bt*sizeof(PetscScalar));CHKERRQ(ierr); 902674ae819SStefano Zampini /* Store upper triangular part of correlation matrix */ 903674ae819SStefano Zampini for (j=0;j<temp_constraints;j++) { 904674ae819SStefano Zampini for (k=0;k<j+1;k++) { 9055b08dc53SStefano Zampini PetscStackCallBLAS("BLASdot",correlation_mat[j*temp_constraints+k]=BLASdot_(&Bs,&temp_quadrature_constraint[temp_indices[temp_start_ptr+j]],&Bone,&temp_quadrature_constraint[temp_indices[temp_start_ptr+k]],&Bone)); 906674ae819SStefano Zampini 907674ae819SStefano Zampini } 908674ae819SStefano Zampini } 909674ae819SStefano Zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 910674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX) 911674ae819SStefano Zampini /* LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,work,&lwork,&lierr); */ 912f7c40c41SStefano Zampini PetscStackCallBLAS("LAPACKsyevx",LAPACKsyevx_("V","A","U",&Bt,correlation_mat,&Bt,&dummy_real,&dummy_real,&dummy_int,&dummy_int,&abs_tol,&eigs_found,singular_vals,singular_vectors,&Bt,work,&lwork,iwork,ifail,&lierr)); 913674ae819SStefano Zampini #else 914674ae819SStefano Zampini /* LAPACK call is missing here! TODO */ 915674ae819SStefano Zampini SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for complexes when PETSC_MISSING_GESVD = 1"); 916674ae819SStefano Zampini #endif 917674ae819SStefano Zampini if (lierr) { 918674ae819SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYEVX Lapack routine %d",(int)lierr); 919674ae819SStefano Zampini } 920674ae819SStefano Zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 921674ae819SStefano Zampini /* retain eigenvalues greater than tol: note that lapack SYEV gives eigs in ascending order */ 922674ae819SStefano Zampini j=0; 923674ae819SStefano Zampini while (j < Bt && singular_vals[j] < tol) j++; 924674ae819SStefano Zampini total_counts=total_counts-j; 925674ae819SStefano Zampini if (j<temp_constraints) { 926674ae819SStefano Zampini for (k=j;k<Bt;k++) { 927674ae819SStefano Zampini singular_vals[k]=1.0/PetscSqrtReal(singular_vals[k]); 928674ae819SStefano Zampini } 929674ae819SStefano Zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 930f7c40c41SStefano Zampini PetscStackCallBLAS("BLASgemm_",BLASgemm_("N","N",&Bs,&Bt,&Bt,&one,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,correlation_mat,&Bt,&zero,temp_basis,&Bs)); 931674ae819SStefano Zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 932674ae819SStefano Zampini /* copy POD basis into used quadrature memory */ 933674ae819SStefano Zampini for (k=0;k<Bt-j;k++) { 934674ae819SStefano Zampini for (ii=0;ii<size_of_constraint;ii++) { 935674ae819SStefano Zampini temp_quadrature_constraint[temp_indices[temp_start_ptr+k]+ii]=singular_vals[Bt-1-k]*temp_basis[(Bt-1-k)*size_of_constraint+ii]; 936674ae819SStefano Zampini } 937674ae819SStefano Zampini } 938674ae819SStefano Zampini } 939674ae819SStefano Zampini 940674ae819SStefano Zampini #else /* on missing GESVD */ 941674ae819SStefano Zampini PetscInt min_n = temp_constraints; 942674ae819SStefano Zampini if (min_n > size_of_constraint) min_n = size_of_constraint; 943674ae819SStefano Zampini dummy_int = Bs; 944674ae819SStefano Zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 945674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX) 946f7c40c41SStefano Zampini PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,singular_vals,&dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,work,&lwork,&lierr)); 947674ae819SStefano Zampini #else 948f7c40c41SStefano Zampini PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,singular_vals,&dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,work,&lwork,rwork,&lierr)); 949674ae819SStefano Zampini #endif 950674ae819SStefano Zampini if (lierr) { 951674ae819SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SVD Lapack routine %d",(int)lierr); 952674ae819SStefano Zampini } 953674ae819SStefano Zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 954674ae819SStefano Zampini /* retain eigenvalues greater than tol: note that lapack SVD gives eigs in descending order */ 955674ae819SStefano Zampini j = 0; 956674ae819SStefano Zampini while (j < min_n && singular_vals[min_n-j-1] < tol) j++; 957674ae819SStefano Zampini total_counts = total_counts-(PetscInt)Bt+(min_n-j); 958674ae819SStefano Zampini #endif 959674ae819SStefano Zampini } 960674ae819SStefano Zampini } 961674ae819SStefano Zampini /* free index sets of faces, edges and vertices */ 962674ae819SStefano Zampini for (i=0;i<n_ISForFaces;i++) { 963674ae819SStefano Zampini ierr = ISDestroy(&ISForFaces[i]);CHKERRQ(ierr); 964674ae819SStefano Zampini } 965674ae819SStefano Zampini ierr = PetscFree(ISForFaces);CHKERRQ(ierr); 966674ae819SStefano Zampini for (i=0;i<n_ISForEdges;i++) { 967674ae819SStefano Zampini ierr = ISDestroy(&ISForEdges[i]);CHKERRQ(ierr); 968674ae819SStefano Zampini } 969674ae819SStefano Zampini ierr = PetscFree(ISForEdges);CHKERRQ(ierr); 970674ae819SStefano Zampini ierr = ISDestroy(&ISForVertices);CHKERRQ(ierr); 971674ae819SStefano Zampini 972674ae819SStefano Zampini /* set quantities in pcbddc data structure */ 973674ae819SStefano Zampini /* n_vertices defines the number of point primal dofs */ 974674ae819SStefano Zampini /* n_constraints defines the number of averages (they can be point primal dofs if change of basis is requested) */ 975674ae819SStefano Zampini local_primal_size = total_counts; 976674ae819SStefano Zampini pcbddc->n_vertices = n_vertices; 977674ae819SStefano Zampini pcbddc->n_constraints = total_counts-n_vertices; 978674ae819SStefano Zampini pcbddc->local_primal_size = local_primal_size; 979674ae819SStefano Zampini 980674ae819SStefano Zampini /* Create constraint matrix */ 981674ae819SStefano Zampini /* The constraint matrix is used to compute the l2g map of primal dofs */ 982674ae819SStefano Zampini /* so we need to set it up properly either with or without change of basis */ 983674ae819SStefano Zampini ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ConstraintMatrix);CHKERRQ(ierr); 984674ae819SStefano Zampini ierr = MatSetType(pcbddc->ConstraintMatrix,impMatType);CHKERRQ(ierr); 985674ae819SStefano Zampini ierr = MatSetSizes(pcbddc->ConstraintMatrix,local_primal_size,pcis->n,local_primal_size,pcis->n);CHKERRQ(ierr); 986674ae819SStefano Zampini /* compute a local numbering of constraints : vertices first then constraints */ 987674ae819SStefano Zampini ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr); 988674ae819SStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array_vector);CHKERRQ(ierr); 989674ae819SStefano Zampini ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&aux_primal_numbering);CHKERRQ(ierr); 990674ae819SStefano Zampini ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&aux_primal_permutation);CHKERRQ(ierr); 991674ae819SStefano Zampini total_counts=0; 992674ae819SStefano Zampini /* find vertices: subdomain corners plus dofs with basis changed */ 993674ae819SStefano Zampini for (i=0;i<local_primal_size;i++) { 994674ae819SStefano Zampini size_of_constraint=temp_indices[i+1]-temp_indices[i]; 995674ae819SStefano Zampini if (change_basis[i] || size_of_constraint == 1) { 996674ae819SStefano Zampini k=0; 997674ae819SStefano Zampini while(k < size_of_constraint && array_vector[temp_indices_to_constraint[temp_indices[i]+size_of_constraint-k-1]] != 0.0) { 998674ae819SStefano Zampini k=k+1; 999674ae819SStefano Zampini } 1000674ae819SStefano Zampini j=temp_indices_to_constraint[temp_indices[i]+size_of_constraint-k-1]; 1001674ae819SStefano Zampini array_vector[j] = 1.0; 1002674ae819SStefano Zampini aux_primal_numbering[total_counts]=j; 1003674ae819SStefano Zampini aux_primal_permutation[total_counts]=total_counts; 1004674ae819SStefano Zampini total_counts++; 1005674ae819SStefano Zampini } 1006674ae819SStefano Zampini } 1007674ae819SStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array_vector);CHKERRQ(ierr); 1008674ae819SStefano Zampini /* permute indices in order to have a sorted set of vertices */ 1009674ae819SStefano Zampini ierr = PetscSortIntWithPermutation(total_counts,aux_primal_numbering,aux_primal_permutation); 1010674ae819SStefano Zampini /* nonzero structure */ 1011674ae819SStefano Zampini ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 1012674ae819SStefano Zampini for (i=0;i<total_counts;i++) { 1013674ae819SStefano Zampini nnz[i]=1; 1014674ae819SStefano Zampini } 1015674ae819SStefano Zampini j=total_counts; 1016674ae819SStefano Zampini for (i=n_vertices;i<local_primal_size;i++) { 1017674ae819SStefano Zampini if (!change_basis[i]) { 1018674ae819SStefano Zampini nnz[j]=temp_indices[i+1]-temp_indices[i]; 1019674ae819SStefano Zampini j++; 1020674ae819SStefano Zampini } 1021674ae819SStefano Zampini } 1022674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(pcbddc->ConstraintMatrix,0,nnz);CHKERRQ(ierr); 1023674ae819SStefano Zampini ierr = PetscFree(nnz);CHKERRQ(ierr); 1024674ae819SStefano Zampini /* set values in constraint matrix */ 1025674ae819SStefano Zampini for (i=0;i<total_counts;i++) { 1026674ae819SStefano Zampini j = aux_primal_permutation[i]; 1027674ae819SStefano Zampini k = aux_primal_numbering[j]; 1028674ae819SStefano Zampini ierr = MatSetValue(pcbddc->ConstraintMatrix,i,k,1.0,INSERT_VALUES);CHKERRQ(ierr); 1029674ae819SStefano Zampini } 1030674ae819SStefano Zampini for (i=n_vertices;i<local_primal_size;i++) { 1031674ae819SStefano Zampini if (!change_basis[i]) { 1032674ae819SStefano Zampini size_of_constraint=temp_indices[i+1]-temp_indices[i]; 1033674ae819SStefano Zampini ierr = MatSetValues(pcbddc->ConstraintMatrix,1,&total_counts,size_of_constraint,&temp_indices_to_constraint[temp_indices[i]],&temp_quadrature_constraint[temp_indices[i]],INSERT_VALUES);CHKERRQ(ierr); 1034674ae819SStefano Zampini total_counts++; 1035674ae819SStefano Zampini } 1036674ae819SStefano Zampini } 1037674ae819SStefano Zampini ierr = PetscFree(aux_primal_numbering);CHKERRQ(ierr); 1038674ae819SStefano Zampini ierr = PetscFree(aux_primal_permutation);CHKERRQ(ierr); 1039674ae819SStefano Zampini /* assembling */ 1040674ae819SStefano Zampini ierr = MatAssemblyBegin(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1041674ae819SStefano Zampini ierr = MatAssemblyEnd(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1042674ae819SStefano Zampini 1043674ae819SStefano Zampini /* Create matrix for change of basis. We don't need it in case pcbddc->use_change_of_basis is FALSE */ 1044674ae819SStefano Zampini if (pcbddc->use_change_of_basis) { 1045674ae819SStefano Zampini ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr); 1046674ae819SStefano Zampini ierr = MatSetType(pcbddc->ChangeOfBasisMatrix,impMatType);CHKERRQ(ierr); 1047674ae819SStefano Zampini ierr = MatSetSizes(pcbddc->ChangeOfBasisMatrix,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B);CHKERRQ(ierr); 1048674ae819SStefano Zampini /* work arrays */ 1049674ae819SStefano Zampini /* we need to reuse these arrays, so we free them */ 1050674ae819SStefano Zampini ierr = PetscFree(temp_basis);CHKERRQ(ierr); 1051674ae819SStefano Zampini ierr = PetscFree(work);CHKERRQ(ierr); 1052674ae819SStefano Zampini ierr = PetscMalloc(pcis->n_B*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 1053674ae819SStefano Zampini ierr = PetscMalloc((nnsp_addone+nnsp_size)*(nnsp_addone+nnsp_size)*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr); 1054674ae819SStefano Zampini ierr = PetscMalloc((nnsp_addone+nnsp_size)*sizeof(PetscScalar),&work);CHKERRQ(ierr); 1055674ae819SStefano Zampini ierr = PetscMalloc((nnsp_addone+nnsp_size)*sizeof(PetscBLASInt),&ipiv);CHKERRQ(ierr); 1056674ae819SStefano Zampini for (i=0;i<pcis->n_B;i++) { 1057674ae819SStefano Zampini nnz[i]=1; 1058674ae819SStefano Zampini } 1059674ae819SStefano Zampini /* Overestimated nonzeros per row */ 1060674ae819SStefano Zampini k=1; 1061674ae819SStefano Zampini for (i=pcbddc->n_vertices;i<local_primal_size;i++) { 1062674ae819SStefano Zampini if (change_basis[i]) { 1063674ae819SStefano Zampini size_of_constraint = temp_indices[i+1]-temp_indices[i]; 1064674ae819SStefano Zampini if (k < size_of_constraint) { 1065674ae819SStefano Zampini k = size_of_constraint; 1066674ae819SStefano Zampini } 1067674ae819SStefano Zampini for (j=0;j<size_of_constraint;j++) { 1068674ae819SStefano Zampini nnz[temp_indices_to_constraint_B[temp_indices[i]+j]] = size_of_constraint; 1069674ae819SStefano Zampini } 1070674ae819SStefano Zampini } 1071674ae819SStefano Zampini } 1072674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(pcbddc->ChangeOfBasisMatrix,0,nnz);CHKERRQ(ierr); 1073674ae819SStefano Zampini ierr = PetscFree(nnz);CHKERRQ(ierr); 1074674ae819SStefano Zampini /* Temporary array to store indices */ 1075674ae819SStefano Zampini ierr = PetscMalloc(k*sizeof(PetscInt),&is_indices);CHKERRQ(ierr); 1076674ae819SStefano Zampini /* Set initial identity in the matrix */ 1077674ae819SStefano Zampini for (i=0;i<pcis->n_B;i++) { 1078674ae819SStefano Zampini ierr = MatSetValue(pcbddc->ChangeOfBasisMatrix,i,i,1.0,INSERT_VALUES);CHKERRQ(ierr); 1079674ae819SStefano Zampini } 1080674ae819SStefano Zampini /* Now we loop on the constraints which need a change of basis */ 1081674ae819SStefano Zampini /* Change of basis matrix is evaluated as the FIRST APPROACH in */ 1082674ae819SStefano Zampini /* Klawonn and Widlund, Dual-primal FETI-DP methods for linear elasticity, (6.2.1) */ 1083674ae819SStefano Zampini temp_constraints = 0; 1084674ae819SStefano Zampini if (pcbddc->n_vertices < local_primal_size) { 1085674ae819SStefano Zampini temp_start_ptr = temp_indices_to_constraint_B[temp_indices[pcbddc->n_vertices]]; 1086674ae819SStefano Zampini } 1087674ae819SStefano Zampini for (i=pcbddc->n_vertices;i<local_primal_size;i++) { 1088674ae819SStefano Zampini if (change_basis[i]) { 1089674ae819SStefano Zampini compute_submatrix = PETSC_FALSE; 1090674ae819SStefano Zampini useksp = PETSC_FALSE; 1091674ae819SStefano Zampini if (temp_start_ptr == temp_indices_to_constraint_B[temp_indices[i]]) { 1092674ae819SStefano Zampini temp_constraints++; 1093674ae819SStefano Zampini if (i == local_primal_size -1 || temp_start_ptr != temp_indices_to_constraint_B[temp_indices[i+1]]) { 1094674ae819SStefano Zampini compute_submatrix = PETSC_TRUE; 1095674ae819SStefano Zampini } 1096674ae819SStefano Zampini } 1097674ae819SStefano Zampini if (compute_submatrix) { 1098674ae819SStefano Zampini if (temp_constraints > 1 || pcbddc->use_nnsp_true) { 1099674ae819SStefano Zampini useksp = PETSC_TRUE; 1100674ae819SStefano Zampini } 1101674ae819SStefano Zampini size_of_constraint = temp_indices[i+1]-temp_indices[i]; 1102674ae819SStefano Zampini if (useksp) { /* experimental TODO: reuse KSP and MAT instead of creating them each time */ 1103674ae819SStefano Zampini ierr = MatCreate(PETSC_COMM_SELF,&temp_mat);CHKERRQ(ierr); 1104674ae819SStefano Zampini ierr = MatSetType(temp_mat,impMatType);CHKERRQ(ierr); 1105674ae819SStefano Zampini ierr = MatSetSizes(temp_mat,size_of_constraint,size_of_constraint,size_of_constraint,size_of_constraint);CHKERRQ(ierr); 1106674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(temp_mat,size_of_constraint,NULL);CHKERRQ(ierr); 1107674ae819SStefano Zampini } 1108674ae819SStefano Zampini /* First _size_of_constraint-temp_constraints_ columns */ 1109674ae819SStefano Zampini dual_dofs = size_of_constraint-temp_constraints; 1110674ae819SStefano Zampini start_constraint = i+1-temp_constraints; 1111674ae819SStefano Zampini for (s=0;s<dual_dofs;s++) { 1112674ae819SStefano Zampini is_indices[0] = s; 1113674ae819SStefano Zampini for (j=0;j<temp_constraints;j++) { 1114674ae819SStefano Zampini for (k=0;k<temp_constraints;k++) { 1115674ae819SStefano Zampini temp_basis[j*temp_constraints+k]=temp_quadrature_constraint[temp_indices[start_constraint+k]+s+j+1]; 1116674ae819SStefano Zampini } 1117674ae819SStefano Zampini work[j]=-temp_quadrature_constraint[temp_indices[start_constraint+j]+s]; 1118674ae819SStefano Zampini is_indices[j+1]=s+j+1; 1119674ae819SStefano Zampini } 1120674ae819SStefano Zampini Bt = temp_constraints; 1121674ae819SStefano Zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1122f7c40c41SStefano Zampini PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&Bt,&Bone,temp_basis,&Bt,ipiv,work,&Bt,&lierr)); 1123674ae819SStefano Zampini if ( lierr ) { 1124674ae819SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GESV Lapack routine %d",(int)lierr); 1125674ae819SStefano Zampini } 1126674ae819SStefano Zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 1127674ae819SStefano Zampini j = temp_indices_to_constraint_B[temp_indices[start_constraint]+s]; 1128674ae819SStefano Zampini ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,temp_constraints,&temp_indices_to_constraint_B[temp_indices[start_constraint]+s+1],1,&j,work,INSERT_VALUES);CHKERRQ(ierr); 1129674ae819SStefano Zampini if (useksp) { 1130674ae819SStefano Zampini /* temp mat with transposed rows and columns */ 1131674ae819SStefano Zampini ierr = MatSetValues(temp_mat,1,&s,temp_constraints,&is_indices[1],work,INSERT_VALUES);CHKERRQ(ierr); 1132674ae819SStefano Zampini ierr = MatSetValue(temp_mat,is_indices[0],is_indices[0],1.0,INSERT_VALUES);CHKERRQ(ierr); 1133674ae819SStefano Zampini } 1134674ae819SStefano Zampini } 1135674ae819SStefano Zampini if (useksp) { 1136674ae819SStefano Zampini /* last rows of temp_mat */ 1137674ae819SStefano Zampini for (j=0;j<size_of_constraint;j++) { 1138674ae819SStefano Zampini is_indices[j] = j; 1139674ae819SStefano Zampini } 1140674ae819SStefano Zampini for (s=0;s<temp_constraints;s++) { 1141674ae819SStefano Zampini k = s + dual_dofs; 1142674ae819SStefano Zampini ierr = MatSetValues(temp_mat,1,&k,size_of_constraint,is_indices,&temp_quadrature_constraint[temp_indices[start_constraint+s]],INSERT_VALUES);CHKERRQ(ierr); 1143674ae819SStefano Zampini } 1144674ae819SStefano Zampini ierr = MatAssemblyBegin(temp_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1145674ae819SStefano Zampini ierr = MatAssemblyEnd(temp_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1146674ae819SStefano Zampini ierr = MatGetVecs(temp_mat,&temp_vec,NULL);CHKERRQ(ierr); 1147674ae819SStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&temp_ksp);CHKERRQ(ierr); 1148674ae819SStefano Zampini ierr = KSPSetOperators(temp_ksp,temp_mat,temp_mat,SAME_PRECONDITIONER);CHKERRQ(ierr); 1149674ae819SStefano Zampini ierr = KSPSetType(temp_ksp,KSPPREONLY);CHKERRQ(ierr); 1150674ae819SStefano Zampini ierr = KSPGetPC(temp_ksp,&temp_pc);CHKERRQ(ierr); 1151674ae819SStefano Zampini ierr = PCSetType(temp_pc,PCLU);CHKERRQ(ierr); 1152674ae819SStefano Zampini ierr = KSPSetUp(temp_ksp);CHKERRQ(ierr); 1153674ae819SStefano Zampini for (s=0;s<temp_constraints;s++) { 1154674ae819SStefano Zampini ierr = VecSet(temp_vec,0.0);CHKERRQ(ierr); 1155674ae819SStefano Zampini ierr = VecSetValue(temp_vec,s+dual_dofs,1.0,INSERT_VALUES);CHKERRQ(ierr); 1156674ae819SStefano Zampini ierr = VecAssemblyBegin(temp_vec);CHKERRQ(ierr); 1157674ae819SStefano Zampini ierr = VecAssemblyEnd(temp_vec);CHKERRQ(ierr); 1158674ae819SStefano Zampini ierr = KSPSolve(temp_ksp,temp_vec,temp_vec);CHKERRQ(ierr); 1159674ae819SStefano Zampini ierr = VecGetArray(temp_vec,&array_vector);CHKERRQ(ierr); 1160674ae819SStefano Zampini j = temp_indices_to_constraint_B[temp_indices[start_constraint+s]+size_of_constraint-s-1]; 1161674ae819SStefano Zampini /* last columns of change of basis matrix associated to new primal dofs */ 1162674ae819SStefano Zampini ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,size_of_constraint,&temp_indices_to_constraint_B[temp_indices[start_constraint+s]],1,&j,array_vector,INSERT_VALUES);CHKERRQ(ierr); 1163674ae819SStefano Zampini ierr = VecRestoreArray(temp_vec,&array_vector);CHKERRQ(ierr); 1164674ae819SStefano Zampini } 1165674ae819SStefano Zampini ierr = MatDestroy(&temp_mat);CHKERRQ(ierr); 1166674ae819SStefano Zampini ierr = KSPDestroy(&temp_ksp);CHKERRQ(ierr); 1167674ae819SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 1168674ae819SStefano Zampini } else { 1169674ae819SStefano Zampini /* last columns of change of basis matrix associated to new primal dofs */ 1170674ae819SStefano Zampini for (s=0;s<temp_constraints;s++) { 1171674ae819SStefano Zampini j = temp_indices_to_constraint_B[temp_indices[start_constraint+s]+size_of_constraint-s-1]; 1172674ae819SStefano Zampini ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,size_of_constraint,&temp_indices_to_constraint_B[temp_indices[start_constraint+s]],1,&j,&temp_quadrature_constraint[temp_indices[start_constraint+s]],INSERT_VALUES);CHKERRQ(ierr); 1173674ae819SStefano Zampini } 1174674ae819SStefano Zampini } 1175674ae819SStefano Zampini /* prepare for the next cycle */ 1176674ae819SStefano Zampini temp_constraints = 0; 1177674ae819SStefano Zampini if (i != local_primal_size -1 ) { 1178674ae819SStefano Zampini temp_start_ptr = temp_indices_to_constraint_B[temp_indices[i+1]]; 1179674ae819SStefano Zampini } 1180674ae819SStefano Zampini } 1181674ae819SStefano Zampini } 1182674ae819SStefano Zampini } 1183674ae819SStefano Zampini /* assembling */ 1184674ae819SStefano Zampini ierr = MatAssemblyBegin(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1185674ae819SStefano Zampini ierr = MatAssemblyEnd(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1186674ae819SStefano Zampini ierr = PetscFree(ipiv);CHKERRQ(ierr); 1187674ae819SStefano Zampini ierr = PetscFree(is_indices);CHKERRQ(ierr); 1188674ae819SStefano Zampini } 1189674ae819SStefano Zampini /* free workspace no longer needed */ 1190674ae819SStefano Zampini ierr = PetscFree(rwork);CHKERRQ(ierr); 1191674ae819SStefano Zampini ierr = PetscFree(work);CHKERRQ(ierr); 1192674ae819SStefano Zampini ierr = PetscFree(temp_basis);CHKERRQ(ierr); 1193674ae819SStefano Zampini ierr = PetscFree(singular_vals);CHKERRQ(ierr); 1194674ae819SStefano Zampini ierr = PetscFree(correlation_mat);CHKERRQ(ierr); 1195674ae819SStefano Zampini ierr = PetscFree(temp_indices);CHKERRQ(ierr); 1196674ae819SStefano Zampini ierr = PetscFree(change_basis);CHKERRQ(ierr); 1197674ae819SStefano Zampini ierr = PetscFree(temp_indices_to_constraint);CHKERRQ(ierr); 1198674ae819SStefano Zampini ierr = PetscFree(temp_indices_to_constraint_B);CHKERRQ(ierr); 1199674ae819SStefano Zampini ierr = PetscFree(local_to_B);CHKERRQ(ierr); 1200674ae819SStefano Zampini ierr = PetscFree(temp_quadrature_constraint);CHKERRQ(ierr); 1201674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD) 1202674ae819SStefano Zampini ierr = PetscFree(iwork);CHKERRQ(ierr); 1203674ae819SStefano Zampini ierr = PetscFree(ifail);CHKERRQ(ierr); 1204674ae819SStefano Zampini ierr = PetscFree(singular_vectors);CHKERRQ(ierr); 1205674ae819SStefano Zampini #endif 1206674ae819SStefano Zampini for (k=0;k<nnsp_size;k++) { 1207674ae819SStefano Zampini ierr = VecDestroy(&localnearnullsp[k]);CHKERRQ(ierr); 1208674ae819SStefano Zampini } 1209674ae819SStefano Zampini ierr = PetscFree(localnearnullsp);CHKERRQ(ierr); 1210674ae819SStefano Zampini PetscFunctionReturn(0); 1211674ae819SStefano Zampini } 1212674ae819SStefano Zampini 1213674ae819SStefano Zampini #undef __FUNCT__ 1214674ae819SStefano Zampini #define __FUNCT__ "PCBDDCAnalyzeInterface" 1215674ae819SStefano Zampini PetscErrorCode PCBDDCAnalyzeInterface(PC pc) 1216674ae819SStefano Zampini { 1217674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1218674ae819SStefano Zampini PC_IS *pcis = (PC_IS*)pc->data; 1219674ae819SStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 1220674ae819SStefano Zampini PetscInt bs,ierr,i,vertex_size; 1221674ae819SStefano Zampini PetscViewer viewer=pcbddc->dbg_viewer; 1222674ae819SStefano Zampini 1223674ae819SStefano Zampini PetscFunctionBegin; 1224674ae819SStefano Zampini /* Init local Graph struct */ 1225674ae819SStefano Zampini ierr = PCBDDCGraphInit(pcbddc->mat_graph,matis->mapping);CHKERRQ(ierr); 1226674ae819SStefano Zampini 1227575ad6abSStefano Zampini /* Check validity of the csr graph passed in by the user */ 1228575ad6abSStefano Zampini if (pcbddc->mat_graph->nvtxs_csr != pcbddc->mat_graph->nvtxs) { 1229575ad6abSStefano Zampini ierr = PCBDDCGraphResetCSR(pcbddc->mat_graph);CHKERRQ(ierr); 1230575ad6abSStefano Zampini } 1231674ae819SStefano Zampini /* Set default CSR adjacency of local dofs if not provided by the user with PCBDDCSetLocalAdjacencyGraph */ 1232674ae819SStefano Zampini if (!pcbddc->mat_graph->xadj || !pcbddc->mat_graph->adjncy) { 1233674ae819SStefano Zampini Mat mat_adj; 1234674ae819SStefano Zampini const PetscInt *xadj,*adjncy; 1235674ae819SStefano Zampini PetscBool flg_row=PETSC_TRUE; 1236674ae819SStefano Zampini 1237674ae819SStefano Zampini ierr = MatConvert(matis->A,MATMPIADJ,MAT_INITIAL_MATRIX,&mat_adj);CHKERRQ(ierr); 1238674ae819SStefano Zampini ierr = MatGetRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr); 1239674ae819SStefano Zampini if (!flg_row) { 1240674ae819SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatGetRowIJ called in %s\n",__FUNCT__); 1241674ae819SStefano Zampini } 1242674ae819SStefano Zampini ierr = PCBDDCSetLocalAdjacencyGraph(pc,i,xadj,adjncy,PETSC_COPY_VALUES);CHKERRQ(ierr); 1243674ae819SStefano Zampini ierr = MatRestoreRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr); 1244674ae819SStefano Zampini if (!flg_row) { 1245674ae819SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 1246674ae819SStefano Zampini } 1247674ae819SStefano Zampini ierr = MatDestroy(&mat_adj);CHKERRQ(ierr); 1248674ae819SStefano Zampini } 1249674ae819SStefano Zampini 1250674ae819SStefano Zampini /* Set default dofs' splitting if no information has been provided by the user with PCBDDCSetDofsSplitting */ 1251674ae819SStefano Zampini vertex_size = 1; 1252674ae819SStefano Zampini if (!pcbddc->n_ISForDofs) { 1253674ae819SStefano Zampini IS *custom_ISForDofs; 1254674ae819SStefano Zampini 1255674ae819SStefano Zampini ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 1256674ae819SStefano Zampini ierr = PetscMalloc(bs*sizeof(IS),&custom_ISForDofs);CHKERRQ(ierr); 1257674ae819SStefano Zampini for (i=0;i<bs;i++) { 1258674ae819SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,pcis->n/bs,i,bs,&custom_ISForDofs[i]);CHKERRQ(ierr); 1259674ae819SStefano Zampini } 1260674ae819SStefano Zampini ierr = PCBDDCSetDofsSplitting(pc,bs,custom_ISForDofs);CHKERRQ(ierr); 1261674ae819SStefano Zampini /* remove my references to IS objects */ 1262674ae819SStefano Zampini for (i=0;i<bs;i++) { 1263674ae819SStefano Zampini ierr = ISDestroy(&custom_ISForDofs[i]);CHKERRQ(ierr); 1264674ae819SStefano Zampini } 1265674ae819SStefano Zampini ierr = PetscFree(custom_ISForDofs);CHKERRQ(ierr); 1266674ae819SStefano Zampini } else { /* mat block size as vertex size (used for elasticity) */ 1267674ae819SStefano Zampini ierr = MatGetBlockSize(matis->A,&vertex_size);CHKERRQ(ierr); 1268674ae819SStefano Zampini } 1269674ae819SStefano Zampini 1270674ae819SStefano Zampini /* Setup of Graph */ 1271674ae819SStefano Zampini ierr = PCBDDCGraphSetUp(pcbddc->mat_graph,vertex_size,pcbddc->NeumannBoundaries,pcbddc->DirichletBoundaries,pcbddc->n_ISForDofs,pcbddc->ISForDofs,pcbddc->user_primal_vertices); 1272674ae819SStefano Zampini 1273674ae819SStefano Zampini /* Graph's connected components analysis */ 1274674ae819SStefano Zampini ierr = PCBDDCGraphComputeConnectedComponents(pcbddc->mat_graph);CHKERRQ(ierr); 1275674ae819SStefano Zampini 1276674ae819SStefano Zampini /* print some info to stdout */ 1277674ae819SStefano Zampini if (pcbddc->dbg_flag) { 1278e49050b4SStefano Zampini ierr = PCBDDCGraphASCIIView(pcbddc->mat_graph,pcbddc->dbg_flag,viewer); 1279674ae819SStefano Zampini } 1280674ae819SStefano Zampini PetscFunctionReturn(0); 1281674ae819SStefano Zampini } 1282674ae819SStefano Zampini 1283674ae819SStefano Zampini #undef __FUNCT__ 1284674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGetPrimalVerticesLocalIdx" 1285674ae819SStefano Zampini PetscErrorCode PCBDDCGetPrimalVerticesLocalIdx(PC pc, PetscInt *n_vertices, PetscInt *vertices_idx[]) 1286674ae819SStefano Zampini { 1287674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 1288674ae819SStefano Zampini PetscInt *vertices,*row_cmat_indices,n,i,size_of_constraint,local_primal_size; 1289674ae819SStefano Zampini PetscErrorCode ierr; 1290674ae819SStefano Zampini 1291674ae819SStefano Zampini PetscFunctionBegin; 1292674ae819SStefano Zampini n = 0; 1293674ae819SStefano Zampini vertices = 0; 1294674ae819SStefano Zampini if (pcbddc->ConstraintMatrix) { 1295674ae819SStefano Zampini ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&i);CHKERRQ(ierr); 1296b120a5c6SStefano Zampini for (i=0;i<local_primal_size;i++) { 1297b120a5c6SStefano Zampini ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 1298b120a5c6SStefano Zampini if (size_of_constraint == 1) n++; 1299b120a5c6SStefano Zampini ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 1300b120a5c6SStefano Zampini } 1301b120a5c6SStefano Zampini ierr = PetscMalloc(n*sizeof(PetscInt),&vertices);CHKERRQ(ierr); 1302b120a5c6SStefano Zampini n = 0; 1303674ae819SStefano Zampini for (i=0;i<local_primal_size;i++) { 1304674ae819SStefano Zampini ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 1305674ae819SStefano Zampini if (size_of_constraint == 1) { 1306674ae819SStefano Zampini vertices[n++]=row_cmat_indices[0]; 1307674ae819SStefano Zampini } 1308674ae819SStefano Zampini ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 1309674ae819SStefano Zampini } 1310674ae819SStefano Zampini } 1311674ae819SStefano Zampini *n_vertices = n; 1312674ae819SStefano Zampini *vertices_idx = vertices; 1313674ae819SStefano Zampini PetscFunctionReturn(0); 1314674ae819SStefano Zampini } 1315674ae819SStefano Zampini 1316674ae819SStefano Zampini #undef __FUNCT__ 1317674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGetPrimalConstraintsLocalIdx" 1318674ae819SStefano Zampini PetscErrorCode PCBDDCGetPrimalConstraintsLocalIdx(PC pc, PetscInt *n_constraints, PetscInt *constraints_idx[]) 1319674ae819SStefano Zampini { 1320674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 1321674ae819SStefano Zampini PetscInt *constraints_index,*row_cmat_indices,*row_cmat_global_indices; 1322674ae819SStefano Zampini PetscInt n,i,j,size_of_constraint,local_primal_size,local_size,max_size_of_constraint,min_index,min_loc; 1323674ae819SStefano Zampini PetscBool *touched; 1324674ae819SStefano Zampini PetscErrorCode ierr; 1325674ae819SStefano Zampini 1326674ae819SStefano Zampini PetscFunctionBegin; 1327674ae819SStefano Zampini n = 0; 1328674ae819SStefano Zampini constraints_index = 0; 1329674ae819SStefano Zampini if (pcbddc->ConstraintMatrix) { 1330674ae819SStefano Zampini ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&local_size);CHKERRQ(ierr); 1331674ae819SStefano Zampini max_size_of_constraint = 0; 1332674ae819SStefano Zampini for (i=0;i<local_primal_size;i++) { 1333674ae819SStefano Zampini ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 1334674ae819SStefano Zampini if (size_of_constraint > 1) { 1335674ae819SStefano Zampini n++; 1336674ae819SStefano Zampini } 1337674ae819SStefano Zampini max_size_of_constraint = PetscMax(size_of_constraint,max_size_of_constraint); 1338674ae819SStefano Zampini ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 1339674ae819SStefano Zampini } 1340674ae819SStefano Zampini ierr = PetscMalloc(n*sizeof(PetscInt),&constraints_index);CHKERRQ(ierr); 1341674ae819SStefano Zampini ierr = PetscMalloc(max_size_of_constraint*sizeof(PetscInt),&row_cmat_global_indices);CHKERRQ(ierr); 1342674ae819SStefano Zampini ierr = PetscMalloc(local_size*sizeof(PetscBool),&touched);CHKERRQ(ierr); 1343674ae819SStefano Zampini ierr = PetscMemzero(touched,local_size*sizeof(PetscBool));CHKERRQ(ierr); 1344674ae819SStefano Zampini n = 0; 1345674ae819SStefano Zampini for (i=0;i<local_primal_size;i++) { 1346674ae819SStefano Zampini ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 1347674ae819SStefano Zampini if (size_of_constraint > 1) { 1348674ae819SStefano Zampini ierr = ISLocalToGlobalMappingApply(pcbddc->mat_graph->l2gmap,size_of_constraint,row_cmat_indices,row_cmat_global_indices);CHKERRQ(ierr); 1349674ae819SStefano Zampini min_index = row_cmat_global_indices[0]; 1350674ae819SStefano Zampini min_loc = 0; 1351674ae819SStefano Zampini for (j=1;j<size_of_constraint;j++) { 1352674ae819SStefano Zampini /* there can be more than one constraint on a single connected component */ 1353674ae819SStefano Zampini if (min_index > row_cmat_global_indices[j] && !touched[row_cmat_indices[j]]) { 1354674ae819SStefano Zampini min_index = row_cmat_global_indices[j]; 1355674ae819SStefano Zampini min_loc = j; 1356674ae819SStefano Zampini } 1357674ae819SStefano Zampini } 1358674ae819SStefano Zampini touched[row_cmat_indices[min_loc]] = PETSC_TRUE; 1359674ae819SStefano Zampini constraints_index[n++] = row_cmat_indices[min_loc]; 1360674ae819SStefano Zampini } 1361674ae819SStefano Zampini ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 1362674ae819SStefano Zampini } 1363674ae819SStefano Zampini } 1364674ae819SStefano Zampini ierr = PetscFree(touched);CHKERRQ(ierr); 1365674ae819SStefano Zampini ierr = PetscFree(row_cmat_global_indices);CHKERRQ(ierr); 1366674ae819SStefano Zampini *n_constraints = n; 1367674ae819SStefano Zampini *constraints_idx = constraints_index; 1368674ae819SStefano Zampini PetscFunctionReturn(0); 1369674ae819SStefano Zampini } 1370674ae819SStefano Zampini 1371674ae819SStefano Zampini /* the next two functions has been adapted from pcis.c */ 1372674ae819SStefano Zampini #undef __FUNCT__ 1373674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplySchur" 1374674ae819SStefano Zampini PetscErrorCode PCBDDCApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 1375674ae819SStefano Zampini { 1376674ae819SStefano Zampini PetscErrorCode ierr; 1377674ae819SStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 1378674ae819SStefano Zampini 1379674ae819SStefano Zampini PetscFunctionBegin; 1380674ae819SStefano Zampini if (!vec2_B) { vec2_B = v; } 1381674ae819SStefano Zampini ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr); 1382674ae819SStefano Zampini ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr); 1383674ae819SStefano Zampini ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr); 1384674ae819SStefano Zampini ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr); 1385674ae819SStefano Zampini ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr); 1386674ae819SStefano Zampini PetscFunctionReturn(0); 1387674ae819SStefano Zampini } 1388674ae819SStefano Zampini 1389674ae819SStefano Zampini #undef __FUNCT__ 1390674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplySchurTranspose" 1391674ae819SStefano Zampini PetscErrorCode PCBDDCApplySchurTranspose(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 1392674ae819SStefano Zampini { 1393674ae819SStefano Zampini PetscErrorCode ierr; 1394674ae819SStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 1395674ae819SStefano Zampini 1396674ae819SStefano Zampini PetscFunctionBegin; 1397674ae819SStefano Zampini if (!vec2_B) { vec2_B = v; } 1398674ae819SStefano Zampini ierr = MatMultTranspose(pcis->A_BB,v,vec1_B);CHKERRQ(ierr); 1399674ae819SStefano Zampini ierr = MatMultTranspose(pcis->A_BI,v,vec1_D);CHKERRQ(ierr); 1400674ae819SStefano Zampini ierr = KSPSolveTranspose(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr); 1401674ae819SStefano Zampini ierr = MatMultTranspose(pcis->A_IB,vec2_D,vec2_B);CHKERRQ(ierr); 1402674ae819SStefano Zampini ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr); 1403674ae819SStefano Zampini PetscFunctionReturn(0); 1404674ae819SStefano Zampini } 1405674ae819SStefano Zampini 1406674ae819SStefano Zampini #undef __FUNCT__ 1407674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSubsetNumbering" 1408674ae819SStefano Zampini PetscErrorCode PCBDDCSubsetNumbering(MPI_Comm comm,ISLocalToGlobalMapping l2gmap, PetscInt n_local_dofs, PetscInt local_dofs[], PetscInt local_dofs_mult[], PetscInt* n_global_subset, PetscInt* global_numbering_subset[]) 1409674ae819SStefano Zampini { 1410674ae819SStefano Zampini Vec local_vec,global_vec; 1411674ae819SStefano Zampini IS seqis,paris; 1412674ae819SStefano Zampini VecScatter scatter_ctx; 1413674ae819SStefano Zampini PetscScalar *array; 1414674ae819SStefano Zampini PetscInt *temp_global_dofs; 1415674ae819SStefano Zampini PetscScalar globalsum; 1416674ae819SStefano Zampini PetscInt i,j,s; 1417674ae819SStefano Zampini PetscInt nlocals,first_index,old_index,max_local; 1418674ae819SStefano Zampini PetscMPIInt rank_prec_comm,size_prec_comm,max_global; 1419674ae819SStefano Zampini PetscMPIInt *dof_sizes,*dof_displs; 1420674ae819SStefano Zampini PetscBool first_found; 1421674ae819SStefano Zampini PetscErrorCode ierr; 1422674ae819SStefano Zampini 1423674ae819SStefano Zampini PetscFunctionBegin; 1424674ae819SStefano Zampini /* mpi buffers */ 1425674ae819SStefano Zampini MPI_Comm_size(comm,&size_prec_comm); 1426674ae819SStefano Zampini MPI_Comm_rank(comm,&rank_prec_comm); 1427674ae819SStefano Zampini j = ( !rank_prec_comm ? size_prec_comm : 0); 1428674ae819SStefano Zampini ierr = PetscMalloc(j*sizeof(*dof_sizes),&dof_sizes);CHKERRQ(ierr); 1429674ae819SStefano Zampini ierr = PetscMalloc(j*sizeof(*dof_displs),&dof_displs);CHKERRQ(ierr); 1430674ae819SStefano Zampini /* get maximum size of subset */ 1431674ae819SStefano Zampini ierr = PetscMalloc(n_local_dofs*sizeof(PetscInt),&temp_global_dofs);CHKERRQ(ierr); 1432674ae819SStefano Zampini ierr = ISLocalToGlobalMappingApply(l2gmap,n_local_dofs,local_dofs,temp_global_dofs);CHKERRQ(ierr); 1433674ae819SStefano Zampini max_local = 0; 1434674ae819SStefano Zampini if (n_local_dofs) { 1435674ae819SStefano Zampini max_local = temp_global_dofs[0]; 1436674ae819SStefano Zampini for (i=1;i<n_local_dofs;i++) { 1437674ae819SStefano Zampini if (max_local < temp_global_dofs[i] ) { 1438674ae819SStefano Zampini max_local = temp_global_dofs[i]; 1439674ae819SStefano Zampini } 1440674ae819SStefano Zampini } 1441674ae819SStefano Zampini } 1442674ae819SStefano Zampini ierr = MPI_Allreduce(&max_local,&max_global,1,MPIU_INT,MPI_MAX,comm); 1443674ae819SStefano Zampini max_global++; 1444674ae819SStefano Zampini max_local = 0; 1445674ae819SStefano Zampini if (n_local_dofs) { 1446674ae819SStefano Zampini max_local = local_dofs[0]; 1447674ae819SStefano Zampini for (i=1;i<n_local_dofs;i++) { 1448674ae819SStefano Zampini if (max_local < local_dofs[i] ) { 1449674ae819SStefano Zampini max_local = local_dofs[i]; 1450674ae819SStefano Zampini } 1451674ae819SStefano Zampini } 1452674ae819SStefano Zampini } 1453674ae819SStefano Zampini max_local++; 1454674ae819SStefano Zampini /* allocate workspace */ 1455674ae819SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr); 1456674ae819SStefano Zampini ierr = VecSetSizes(local_vec,PETSC_DECIDE,max_local);CHKERRQ(ierr); 1457674ae819SStefano Zampini ierr = VecSetType(local_vec,VECSEQ);CHKERRQ(ierr); 1458674ae819SStefano Zampini ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr); 1459674ae819SStefano Zampini ierr = VecSetSizes(global_vec,PETSC_DECIDE,max_global);CHKERRQ(ierr); 1460674ae819SStefano Zampini ierr = VecSetType(global_vec,VECMPI);CHKERRQ(ierr); 1461674ae819SStefano Zampini /* create scatter */ 1462674ae819SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,n_local_dofs,local_dofs,PETSC_COPY_VALUES,&seqis);CHKERRQ(ierr); 1463674ae819SStefano Zampini ierr = ISCreateGeneral(comm,n_local_dofs,temp_global_dofs,PETSC_COPY_VALUES,&paris);CHKERRQ(ierr); 1464674ae819SStefano Zampini ierr = VecScatterCreate(local_vec,seqis,global_vec,paris,&scatter_ctx);CHKERRQ(ierr); 1465674ae819SStefano Zampini ierr = ISDestroy(&seqis);CHKERRQ(ierr); 1466674ae819SStefano Zampini ierr = ISDestroy(&paris);CHKERRQ(ierr); 1467674ae819SStefano Zampini /* init array */ 1468674ae819SStefano Zampini ierr = VecSet(global_vec,0.0);CHKERRQ(ierr); 1469674ae819SStefano Zampini ierr = VecSet(local_vec,0.0);CHKERRQ(ierr); 1470674ae819SStefano Zampini ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 1471674ae819SStefano Zampini if (local_dofs_mult) { 1472674ae819SStefano Zampini for (i=0;i<n_local_dofs;i++) { 1473674ae819SStefano Zampini array[local_dofs[i]]=(PetscScalar)local_dofs_mult[i]; 1474674ae819SStefano Zampini } 1475674ae819SStefano Zampini } else { 1476674ae819SStefano Zampini for (i=0;i<n_local_dofs;i++) { 1477674ae819SStefano Zampini array[local_dofs[i]]=1.0; 1478674ae819SStefano Zampini } 1479674ae819SStefano Zampini } 1480674ae819SStefano Zampini ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 1481674ae819SStefano Zampini /* scatter into global vec and get total number of global dofs */ 1482674ae819SStefano Zampini ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1483674ae819SStefano Zampini ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1484674ae819SStefano Zampini ierr = VecSum(global_vec,&globalsum);CHKERRQ(ierr); 14855b08dc53SStefano Zampini *n_global_subset = (PetscInt)PetscRealPart(globalsum); 1486674ae819SStefano Zampini /* Fill global_vec with cumulative function for global numbering */ 1487674ae819SStefano Zampini ierr = VecGetArray(global_vec,&array);CHKERRQ(ierr); 1488674ae819SStefano Zampini ierr = VecGetLocalSize(global_vec,&s);CHKERRQ(ierr); 1489674ae819SStefano Zampini nlocals = 0; 1490674ae819SStefano Zampini first_index = -1; 1491674ae819SStefano Zampini first_found = PETSC_FALSE; 1492674ae819SStefano Zampini for (i=0;i<s;i++) { 14935b08dc53SStefano Zampini if (!first_found && PetscRealPart(array[i]) > 0.0) { 1494674ae819SStefano Zampini first_found = PETSC_TRUE; 1495674ae819SStefano Zampini first_index = i; 1496674ae819SStefano Zampini } 14975b08dc53SStefano Zampini nlocals += (PetscInt)PetscRealPart(array[i]); 1498674ae819SStefano Zampini } 1499674ae819SStefano Zampini ierr = MPI_Gather(&nlocals,1,MPIU_INT,dof_sizes,1,MPIU_INT,0,comm);CHKERRQ(ierr); 1500674ae819SStefano Zampini if (!rank_prec_comm) { 1501674ae819SStefano Zampini dof_displs[0]=0; 1502674ae819SStefano Zampini for (i=1;i<size_prec_comm;i++) { 1503674ae819SStefano Zampini dof_displs[i] = dof_displs[i-1]+dof_sizes[i-1]; 1504674ae819SStefano Zampini } 1505674ae819SStefano Zampini } 1506674ae819SStefano Zampini ierr = MPI_Scatter(dof_displs,1,MPIU_INT,&nlocals,1,MPIU_INT,0,comm);CHKERRQ(ierr); 1507674ae819SStefano Zampini if (first_found) { 1508674ae819SStefano Zampini array[first_index] += (PetscScalar)nlocals; 1509674ae819SStefano Zampini old_index = first_index; 1510674ae819SStefano Zampini for (i=first_index+1;i<s;i++) { 15115b08dc53SStefano Zampini if (PetscRealPart(array[i]) > 0.0) { 1512674ae819SStefano Zampini array[i] += array[old_index]; 1513674ae819SStefano Zampini old_index = i; 1514674ae819SStefano Zampini } 1515674ae819SStefano Zampini } 1516674ae819SStefano Zampini } 1517674ae819SStefano Zampini ierr = VecRestoreArray(global_vec,&array);CHKERRQ(ierr); 1518674ae819SStefano Zampini ierr = VecSet(local_vec,0.0);CHKERRQ(ierr); 1519674ae819SStefano Zampini ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1520674ae819SStefano Zampini ierr = VecScatterEnd (scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1521674ae819SStefano Zampini /* get global ordering of local dofs */ 1522674ae819SStefano Zampini ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 1523674ae819SStefano Zampini if (local_dofs_mult) { 1524674ae819SStefano Zampini for (i=0;i<n_local_dofs;i++) { 15255b08dc53SStefano Zampini temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-local_dofs_mult[i]; 1526674ae819SStefano Zampini } 1527674ae819SStefano Zampini } else { 1528674ae819SStefano Zampini for (i=0;i<n_local_dofs;i++) { 15295b08dc53SStefano Zampini temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-1; 1530674ae819SStefano Zampini } 1531674ae819SStefano Zampini } 1532674ae819SStefano Zampini ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 1533674ae819SStefano Zampini /* free workspace */ 1534674ae819SStefano Zampini ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr); 1535674ae819SStefano Zampini ierr = VecDestroy(&local_vec);CHKERRQ(ierr); 1536674ae819SStefano Zampini ierr = VecDestroy(&global_vec);CHKERRQ(ierr); 1537674ae819SStefano Zampini ierr = PetscFree(dof_sizes);CHKERRQ(ierr); 1538674ae819SStefano Zampini ierr = PetscFree(dof_displs);CHKERRQ(ierr); 1539674ae819SStefano Zampini /* return pointer to global ordering of local dofs */ 1540674ae819SStefano Zampini *global_numbering_subset = temp_global_dofs; 1541674ae819SStefano Zampini PetscFunctionReturn(0); 1542674ae819SStefano Zampini } 1543