Skip to content

Commit a99ed84

Browse files
fixed a problem with phase-change BCasting when option was False
1 parent 48920e6 commit a99ed84

1 file changed

Lines changed: 87 additions & 72 deletions

File tree

src/reader.cpp

Lines changed: 87 additions & 72 deletions
Original file line numberDiff line numberDiff line change
@@ -752,6 +752,18 @@ PetscErrorCode reader(int rank, const char fName[]){
752752
printf("\n V: ");
753753
for (PetscInt i=0;i<n_interfaces+1;i++)
754754
printf("%.3e ",inter_V[i]);
755+
printf("\ncmax: ");
756+
for (PetscInt i=0;i<n_interfaces+1;i++)
757+
printf("%.3e ",cohesion_max[i]);
758+
printf("\ncmin: ");
759+
for (PetscInt i=0;i<n_interfaces+1;i++)
760+
printf("%.3e ",cohesion_min[i]);
761+
printf("\nfmax: ");
762+
for (PetscInt i=0;i<n_interfaces+1;i++)
763+
printf("%.3e ",friction_angle_max[i]);
764+
printf("\nfmin: ");
765+
for (PetscInt i=0;i<n_interfaces+1;i++)
766+
printf("%.3e ",friction_angle_min[i]);
755767
printf("\n\n");
756768

757769
fclose(f_interfaces); // Close file
@@ -775,89 +787,92 @@ PetscErrorCode reader(int rank, const char fName[]){
775787
MPI_Bcast(cohesion_min,n_interfaces+1,MPIU_SCALAR,0,PETSC_COMM_WORLD);
776788
MPI_Bcast(friction_angle_max,n_interfaces+1,MPIU_SCALAR,0,PETSC_COMM_WORLD);
777789
MPI_Bcast(friction_angle_min,n_interfaces+1,MPIU_SCALAR,0,PETSC_COMM_WORLD);
778-
790+
779791
// Read phase change files
780-
PetscCalloc1(n_interfaces+1,&phase_change_unit_number);
781-
for (PetscInt k=0; k<n_interfaces+1; k+=1) phase_change_unit_number[k] = -1;
782-
783-
if ((phase_change==1) && (rank==0))
792+
if (phase_change==1)
784793
{
785-
for (PetscInt i=0; i<n_interfaces+1; i+=1)
794+
PetscCalloc1(n_interfaces+1,&phase_change_unit_number);
795+
for (PetscInt k=0; k<n_interfaces+1; k+=1) phase_change_unit_number[k] = -1;
796+
797+
if (rank==0)
786798
{
787-
char fname[100] = "models/phase_change_";
788-
char str_i[11];
789-
char str_e[10] = ".txt";
790-
791-
sprintf(str_i, "%d", i);
792-
strcat(fname, str_i);
793-
strcat(fname, str_e);
794-
read_phase_change_file(fname, i);
799+
for (PetscInt i=0; i<n_interfaces+1; i+=1)
800+
{
801+
char fname[100] = "models/phase_change_";
802+
char str_i[11];
803+
char str_e[10] = ".txt";
804+
805+
sprintf(str_i, "%d", i);
806+
strcat(fname, str_i);
807+
strcat(fname, str_e);
808+
read_phase_change_file(fname, i);
809+
}
810+
PetscPrintf(PETSC_COMM_WORLD, "\n");
795811
}
796-
PetscPrintf(PETSC_COMM_WORLD, "\n");
797-
}
798812

799-
// Broadcast number of phase files
800-
MPI_Bcast(&num_phase_change_files,1,MPIU_INT,0,PETSC_COMM_WORLD);
813+
// Broadcast number of phase files
814+
MPI_Bcast(&num_phase_change_files,1,MPIU_INT,0,PETSC_COMM_WORLD);
801815

802-
// Broadcast flag files
803-
MPI_Bcast(phase_change_unit_number,n_interfaces+1,MPIU_INT,0,PETSC_COMM_WORLD);
816+
// Broadcast flag files
817+
MPI_Bcast(phase_change_unit_number,n_interfaces+1,MPIU_INT,0,PETSC_COMM_WORLD);
804818

805-
// Allocate size arrays
806-
PetscCalloc1(num_phase_change_files+1,&p_size);
807-
PetscCalloc1(num_phase_change_files+1,&p_cum_size);
808-
PetscCalloc1(num_phase_change_files+1,&t_size);
809-
PetscCalloc1(num_phase_change_files+1,&t_cum_size);
810-
PetscCalloc1(num_phase_change_files+1,&d_size);
811-
PetscCalloc1(num_phase_change_files+1,&d_cum_size);
819+
// Allocate size arrays
820+
PetscCalloc1(num_phase_change_files+1,&p_size);
821+
PetscCalloc1(num_phase_change_files+1,&p_cum_size);
822+
PetscCalloc1(num_phase_change_files+1,&t_size);
823+
PetscCalloc1(num_phase_change_files+1,&t_cum_size);
824+
PetscCalloc1(num_phase_change_files+1,&d_size);
825+
PetscCalloc1(num_phase_change_files+1,&d_cum_size);
812826

813-
// Copy from rank 0 to global arrays
814-
if (rank==0)
815-
{
816-
memcpy(p_size,p_size_0,(num_phase_change_files+1)*sizeof(PetscScalar));
817-
memcpy(p_cum_size,p_cum_size_0,(num_phase_change_files+1)*sizeof(PetscScalar));
818-
memcpy(t_size,t_size_0,(num_phase_change_files+1)*sizeof(PetscScalar));
819-
memcpy(t_cum_size,t_cum_size_0,(num_phase_change_files+1)*sizeof(PetscScalar));
820-
memcpy(d_size,d_size_0,(num_phase_change_files+1)*sizeof(PetscScalar));
821-
memcpy(d_cum_size,d_cum_size_0,(num_phase_change_files+1)*sizeof(PetscScalar));
822-
}
827+
// Copy from rank 0 to global arrays
828+
if (rank==0)
829+
{
830+
memcpy(p_size,p_size_0,(num_phase_change_files+1)*sizeof(PetscScalar));
831+
memcpy(p_cum_size,p_cum_size_0,(num_phase_change_files+1)*sizeof(PetscScalar));
832+
memcpy(t_size,t_size_0,(num_phase_change_files+1)*sizeof(PetscScalar));
833+
memcpy(t_cum_size,t_cum_size_0,(num_phase_change_files+1)*sizeof(PetscScalar));
834+
memcpy(d_size,d_size_0,(num_phase_change_files+1)*sizeof(PetscScalar));
835+
memcpy(d_cum_size,d_cum_size_0,(num_phase_change_files+1)*sizeof(PetscScalar));
836+
}
823837

824-
// Broadcast size arrays
825-
MPI_Bcast(p_size,num_phase_change_files+1,MPIU_INT,0,PETSC_COMM_WORLD);
826-
MPI_Bcast(p_cum_size,num_phase_change_files+1,MPIU_INT,0,PETSC_COMM_WORLD);
827-
MPI_Bcast(t_size,num_phase_change_files+1,MPIU_INT,0,PETSC_COMM_WORLD);
828-
MPI_Bcast(t_cum_size,num_phase_change_files+1,MPIU_INT,0,PETSC_COMM_WORLD);
829-
MPI_Bcast(d_size,num_phase_change_files+1,MPIU_INT,0,PETSC_COMM_WORLD);
830-
MPI_Bcast(d_cum_size,num_phase_change_files+1,MPIU_INT,0,PETSC_COMM_WORLD);
831-
832-
// Allocate data arrays
833-
PetscCalloc1(p_cum_size[num_phase_change_files],&phase_pressure);
834-
PetscCalloc1(t_cum_size[num_phase_change_files],&phase_temperature);
835-
PetscCalloc1(d_cum_size[num_phase_change_files],&phase_density);
836-
837-
// Copy from rank 0 to global data arrays
838-
if (rank==0)
839-
{
840-
memcpy(phase_pressure,phase_pressure_0,(p_cum_size[num_phase_change_files])*sizeof(PetscScalar));
841-
memcpy(phase_temperature,phase_temperature_0,(t_cum_size[num_phase_change_files])*sizeof(PetscScalar));
842-
memcpy(phase_density,phase_density_0,(d_cum_size[num_phase_change_files])*sizeof(PetscScalar));
843-
}
838+
// Broadcast size arrays
839+
MPI_Bcast(p_size,num_phase_change_files+1,MPIU_INT,0,PETSC_COMM_WORLD);
840+
MPI_Bcast(p_cum_size,num_phase_change_files+1,MPIU_INT,0,PETSC_COMM_WORLD);
841+
MPI_Bcast(t_size,num_phase_change_files+1,MPIU_INT,0,PETSC_COMM_WORLD);
842+
MPI_Bcast(t_cum_size,num_phase_change_files+1,MPIU_INT,0,PETSC_COMM_WORLD);
843+
MPI_Bcast(d_size,num_phase_change_files+1,MPIU_INT,0,PETSC_COMM_WORLD);
844+
MPI_Bcast(d_cum_size,num_phase_change_files+1,MPIU_INT,0,PETSC_COMM_WORLD);
845+
846+
// Allocate data arrays
847+
PetscCalloc1(p_cum_size[num_phase_change_files],&phase_pressure);
848+
PetscCalloc1(t_cum_size[num_phase_change_files],&phase_temperature);
849+
PetscCalloc1(d_cum_size[num_phase_change_files],&phase_density);
850+
851+
// Copy from rank 0 to global data arrays
852+
if (rank==0)
853+
{
854+
memcpy(phase_pressure,phase_pressure_0,(p_cum_size[num_phase_change_files])*sizeof(PetscScalar));
855+
memcpy(phase_temperature,phase_temperature_0,(t_cum_size[num_phase_change_files])*sizeof(PetscScalar));
856+
memcpy(phase_density,phase_density_0,(d_cum_size[num_phase_change_files])*sizeof(PetscScalar));
857+
}
844858

845-
// Broadcast data arrays
846-
MPI_Bcast(phase_pressure,p_cum_size[num_phase_change_files],MPIU_SCALAR,0,PETSC_COMM_WORLD);
847-
MPI_Bcast(phase_temperature,t_cum_size[num_phase_change_files],MPIU_SCALAR,0,PETSC_COMM_WORLD);
848-
MPI_Bcast(phase_density,d_cum_size[num_phase_change_files],MPIU_SCALAR,0,PETSC_COMM_WORLD);
849-
850-
// Free rank 0
851-
PetscFree(p_size_0);
852-
PetscFree(p_cum_size_0);
853-
PetscFree(t_size_0);
854-
PetscFree(t_cum_size_0);
855-
PetscFree(d_size_0);
856-
PetscFree(d_cum_size_0);
857-
PetscFree(phase_pressure_0);
858-
PetscFree(phase_temperature_0);
859-
PetscFree(phase_density_0);
859+
// Broadcast data arrays
860+
MPI_Bcast(phase_pressure,p_cum_size[num_phase_change_files],MPIU_SCALAR,0,PETSC_COMM_WORLD);
861+
MPI_Bcast(phase_temperature,t_cum_size[num_phase_change_files],MPIU_SCALAR,0,PETSC_COMM_WORLD);
862+
MPI_Bcast(phase_density,d_cum_size[num_phase_change_files],MPIU_SCALAR,0,PETSC_COMM_WORLD);
863+
864+
// Free rank 0
865+
PetscFree(p_size_0);
866+
PetscFree(p_cum_size_0);
867+
PetscFree(t_size_0);
868+
PetscFree(t_cum_size_0);
869+
PetscFree(d_size_0);
870+
PetscFree(d_cum_size_0);
871+
PetscFree(phase_pressure_0);
872+
PetscFree(phase_temperature_0);
873+
PetscFree(phase_density_0);
860874

875+
}
861876
// Broadcast, special cases
862877
// MPI_Bcast(&n_interfaces,1,MPI_INT,0,PETSC_COMM_WORLD); // Broadcast after interfaces.txt
863878

0 commit comments

Comments
 (0)